Public API
QInchworm
— ModuleA quasi Monte Carlo inchworm impurity solver for multi-orbital fermionic models.
QInchworm.expansion
QInchworm.expansion
— ModuleStrong coupling pseudo-particle expansion problem.
Exports
QInchworm.expansion.Expansion
— Typestruct Expansion{ScalarGF<:Keldysh.AbstractTimeGF{ComplexF64, true}, PPGF_P0<:Union{Vector{Keldysh.GenericTimeGF{ComplexF64, false, Keldysh.FullTimeGrid}}, Vector{Keldysh.ImaginaryTimeGF{ComplexF64, false}}, Vector{QInchworm.exact_atomic_ppgf.ExactAtomicPPGF}, Vector{QInchworm.spline_gf.IncSplineImaginaryTimeGF{ComplexF64, false}}, Vector{QInchworm.spline_gf.SplineInterpolatedGF{Keldysh.GenericTimeGF{ComplexF64, false, Keldysh.FullTimeGrid}, ComplexF64, false}}, Vector{QInchworm.spline_gf.SplineInterpolatedGF{Keldysh.ImaginaryTimeGF{ComplexF64, false}, ComplexF64, false}}}, PPGF_P<:Union{Vector{Keldysh.GenericTimeGF{ComplexF64, false, Keldysh.FullTimeGrid}}, Vector{Keldysh.ImaginaryTimeGF{ComplexF64, false}}, Vector{QInchworm.exact_atomic_ppgf.ExactAtomicPPGF}, Vector{QInchworm.spline_gf.IncSplineImaginaryTimeGF{ComplexF64, false}}, Vector{QInchworm.spline_gf.SplineInterpolatedGF{Keldysh.GenericTimeGF{ComplexF64, false, Keldysh.FullTimeGrid}, ComplexF64, false}}, Vector{QInchworm.spline_gf.SplineInterpolatedGF{Keldysh.ImaginaryTimeGF{ComplexF64, false}, ComplexF64, false}}}}
The Expansion
structure contains the components needed to define a strong coupling pseudo-particle expansion problem.
Fields
ed
: Exact diagonalization solver for the local degrees of freedomP0
: Non-interacting propagator (pseudo-particle Green's function)P
: Interacting propagator (pseudo-particle Green's function)pairs
: List of pseudo-particle interactionsdeterminants
: List of hybridization function determinants (not implemented yet)corr_operators
: List of operator pairs used in accumulation of two-point correlation functionsidentity_mat
: Block matrix representation of the identity operatorpair_operator_mat
: Block matrix representation of paired operators (operator_i
,operator_f
)corr_operators_mat
: Block matrix representation ofcorr_operators
subspace_attachable_pairs
: Given a subspace indexs_i
, lists indices of all interaction pairs with operator_i acting non-trivially on that subspace.
QInchworm.expansion.Expansion
— MethodExpansion(
ed::KeldyshED.EDCore,
grid::Keldysh.AbstractTimeGrid,
interaction_pairs::Array{QInchworm.expansion.InteractionPair{ScalarGF}, 1};
corr_operators,
interpolate_ppgf
) -> QInchworm.expansion.Expansion
Parameters
ed
: Exact diagonalization solution of the local problem.grid
: Contour time grid to define the local propagators on.interaction_pairs
: The list of pair interactions to expand in.corr_operators
: The list of operator pairs used in accumulation of two-point correlation functions.interpolate_ppgf
: Use a quadratic spline interpolation to represent and evaluate the local propagators. Currently works only with the imaginary time propagators.
QInchworm.expansion.Expansion
— MethodExpansion(
hamiltonian::KeldyshED.Operators.OperatorExpr,
soi::KeldyshED.Hilbert.SetOfIndices,
grid::Keldysh.ImaginaryTimeGrid;
hybridization,
nn_interaction,
corr_operators,
interpolate_ppgf
)
A higher-level constructor of Expansion
that solves the local problem defined by a Hamiltonian and internally generates a list of pseudo-particle pair interactions from hybridization and $nn$-interaction functions.
Parameters
hamiltonian
: Hamiltonian of the local problem.soi
: An ordered set of indices carried by creation/annihilation operators of the local problem.grid
: Imaginary time grid to define the local propagators on.hybridization
: A matrix-valued hybridization function $\Delta_{ij}(\tau)$. A correspondence between the matrix elements $(i, j)$ and operators $c^\dagger_i, c_j$ is established bysoi
.nn_interaction
: A matrix-valued $nn$-interaction function $U_{ij}(\tau)$. A correspondence between the matrix elements $(i, j)$ and operators $n_i, n_j$ is established bysoi
.corr_operators
: The list of operator pairs used in accumulation of two-point correlation functions.interpolate_ppgf
: Use a quadratic spline interpolation to represent and evaluate the local propagators.
QInchworm.expansion.InteractionPair
— Typestruct InteractionPair{ScalarGF<:Keldysh.AbstractTimeGF{ComplexF64, true}}
Data type for pseudo-particle interactions containing two operators and one scalar propagator.
Indexed access to the operators stored in a pair::InteractionPair
is supported: pair[1]
and pair[2]
are equivalent to pair.operator_i
and pair.operator_f
respectively.
Fields
operator_f::KeldyshED.Operators.RealOperatorExpr
: Final time operatoroperator_i::KeldyshED.Operators.RealOperatorExpr
: Initial time operatorpropagator::Keldysh.AbstractTimeGF{ComplexF64, true}
: Scalar propagator
QInchworm.expansion.add_corr_operators!
— Functionadd_corr_operators!(
expansion::QInchworm.expansion.Expansion,
ops::Tuple{KeldyshED.Operators.RealOperatorExpr, KeldyshED.Operators.RealOperatorExpr}
)
Add a pair of operators $(A, B)$ used to measure the two-point correlator $\langle A(t_1) B(t_2)\rangle$ to expansion
.
Parameters
expansion
: Pseudo-particle expansion.ops
: The pair of operators $(A, B)$.
QInchworm.inchworm
QInchworm.inchworm
— ModuleHigh level functions implementing the quasi Monte Carlo inchworm algorithm.
Exports
QInchworm.inchworm.inchworm!
— Functioninchworm!(
expansion::QInchworm.expansion.Expansion,
grid::Keldysh.ImaginaryTimeGrid,
orders,
orders_bare,
N_samples::Int64;
n_pts_after_max,
rand_params,
seq_type
) -> Tuple{Dict, Dict}
Perform a complete qMC inchworm calculation of the bold propagators on the imaginary time segment. Results of the calculation are written into expansion.P
.
Parameters
expansion
: Strong coupling expansion problem.grid
: Imaginary time grid of the bold propagators.orders
: List of expansion orders to be accounted for during a regular inchworm step.orders_bare
: List of expansion orders to be accounted for during the initial inchworm step.N_samples
: Number of samples to be used in qMC integration. Must be a power of 2.n_pts_after_max
: Maximum number of points in the after-$\tau_w$ region to be taken into account. By default, diagrams with all valid numbers of the after-$\tau_w$ points are considered.rand_params
: Parameters of the randomized qMC integration.seq_type
: Type of the (quasi-)random sequence to be used for integration.
Returns
- Order-resolved contributions to the bold propagator as a dictionary
Dict{Int, PPGF}
. - Estimated standard deviations of the order-resolved contributions as a dictionary
Dict{Int, PPGF}
.
QInchworm.inchworm.diff_inchworm!
— Functiondiff_inchworm!(
expansion::QInchworm.expansion.Expansion,
grid::Keldysh.ImaginaryTimeGrid,
orders,
N_samples::Int64;
rand_params,
seq_type
) -> Tuple{Dict, Dict}
Perform a complete qMC inchworm calculation of the bold propagators on the imaginary time segment using the differential formulation of the method described in
"Inchworm Monte Carlo Method for Open Quantum Systems"
Z. Cai, J. Lu and S. Yang
Comm. Pure Appl. Math., 73: 2430-2472 (2020)
Results of the calculation are written into expansion.P
.
Parameters
expansion
: Strong coupling expansion problem.grid
: Imaginary time grid of the bold propagators.orders
: List of expansion orders to be accounted for.N_samples
: Number of samples to be used in qMC integration. Must be a power of 2.rand_params
: Parameters of the randomized qMC integration.seq_type
: Type of the (quasi-)random sequence to be used for integration.
Returns
- Order-resolved contributions to the pseudo-particle self-energy as a dictionary
Dict{Int, PPGF}
. - Estimated standard deviations of the order-resolved contributions as a dictionary
Dict{Int, PPGF}
.
QInchworm.inchworm.correlator_2p
— Functioncorrelator_2p(
expansion::QInchworm.expansion.Expansion,
grid::Keldysh.ImaginaryTimeGrid,
A_B_pair_idx::Int64,
τ::Keldysh.TimeGridPoint,
top_data::Vector{QInchworm.inchworm.TopologiesInputData};
seq_type,
tmr
) -> Tuple{ComplexF64, ComplexF64}
Calculate value of a two-point correlator $\langle A(\tau) B(0)\rangle$ for one value of the imaginary time argument $\tau$. The pair of operators $(A, B)$ used in the calculation is taken from expansion.corr_operators[A_B_pair_idx]
.
Parameters
expansion
: Strong coupling expansion problem.expansion.P
must contain precomputed bold propagators.grid
: Imaginary time grid of the correlator to be computed.A_B_pair_idx
: Index of the $(A, B)$ pair withinexpansion.corr_operators
.τ
: The imaginary time argument $\tau$.top_data
: Accumulation input data.seq_type
: Type of the (quasi-)random sequence to be used for integration.tmr
: ATimerOutput
object used for profiling.
Returns
- Accumulated value of the two-point correlator.
- Estimated standard deviations of the computed correlator.
correlator_2p(
expansion::QInchworm.expansion.Expansion,
grid::Keldysh.ImaginaryTimeGrid,
orders,
N_samples::Int64,
::QInchworm.randomization.RequestStdDev;
rand_params,
seq_type
) -> Tuple{Any, Any}
Calculate a two-point correlator $\langle A(\tau) B(0)\rangle$ on the imaginary time segment. Accumulation is performed for each pair of operators $(A, B)$ in expansion.corr_operators
. Only the operators that are a single monomial in $c/c^\dagger$ are supported.
This method is selected by the flag argument of type RequestStdDev
and returns randomized qMC estimates of both mean and standard deviation of the correlators.
Parameters
expansion
: Strong coupling expansion problem.expansion.P
must contain precomputed bold propagators.grid
: Imaginary time grid of the correlator to be computed.orders
: List of expansion orders to be accounted for.N_samples
: Number of samples to be used in qMC integration. Must be a power of 2.rand_params
: Parameters of the randomized qMC integration.seq_type
: Type of the (quasi-)random sequence to be used for integration.
Returns
- A list of scalar-valued GF objects containing the computed correlators, one element per a pair in
expansion.corr_operators
. - A list of scalar-valued GF objects containing estimated standard deviations of the computed correlators, one element per a pair in
expansion.corr_operators
.
correlator_2p(
expansion::QInchworm.expansion.Expansion,
grid::Keldysh.ImaginaryTimeGrid,
orders,
N_samples::Int64;
rand_params,
seq_type
) -> Any
Calculate a two-point correlator $\langle A(\tau) B(0)\rangle$ on the imaginary time segment. Accumulation is performed for each pair of operators $(A, B)$ in expansion.corr_operators
. Only the operators that are a single monomial in $c/c^\dagger$ are supported.
Parameters
expansion
: Strong coupling expansion problem.expansion.P
must contain precomputed bold propagators.grid
: Imaginary time grid of the correlator to be computed.orders
: List of expansion orders to be accounted for.N_samples
: Number of samples to be used in qMC integration. Must be a power of 2.rand_params
: Parameters of the randomized qMC integration.seq_type
: Type of the (quasi-)random sequence to be used for integration.
Returns
A list of scalar-valued GF objects containing the computed correlators, one element per a pair in expansion.corr_operators
.
QInchworm.randomization
QInchworm.randomization
— ModuleQInchworm.randomization.RandomizationParams
— Typestruct RandomizationParams
Parameters of the randomized qMC integration.
Fields
rng::Union{Nothing, Random.AbstractRNG}
: Random Number Generator used to scramble Sobol sequences ornothing
to disable scrambling
N_seqs::Int64
: Maximal number of scrambled Sobol sequences to be usedtarget_std::Float64
: Target standard deviation of the accumulated quantity (∞-norm for matrices)
QInchworm.randomization.RequestStdDev
— Typestruct RequestStdDev
Singleton type used to select methods that return an estimate of the standard deviation in addition to the mean.
QInchworm.ppgf
QInchworm.ppgf
— ModulePseudo-particle Green's functions (propagators) of finite fermionic systems and related tools.
For a system defined by a time-independent Hamiltonian $\hat H$, the pseudo-particle Green's function (PPGF) is
\[P(z, z') = \left\{ \begin{array}{ll} -i (-1)^{\hat N} e^{-i \hat H(z-z')},& z \succ -i\beta \cap -i\beta \succeq z',\\ -i e^{-i \hat H(z-z')},& \text{otherwise}. \end{array} \right.\]
In particular, on the imaginary time segment alone one has $P(\tau) = -i e^{-\hat H \tau}$.
This operator has a block-diagonal structure determined by the symmetry sectors of $\hat H$, and is stored as a vector of GF containers corresponding to the individual diagonal blocks (FullTimePPGF
, ImaginaryTimePPGF
).
Exports
QInchworm.ppgf.FullTimePPGF
— TypeAn atomic propagator (PPGF) defined on a full Keldysh contour
QInchworm.ppgf.ImaginaryTimePPGF
— TypeAn atomic propagator (PPGF) defined on an imaginary time segment
QInchworm.ppgf.atomic_ppgf
— Functionatomic_ppgf(
β::Float64,
ed::KeldyshED.EDCore
) -> Vector{QInchworm.exact_atomic_ppgf.ExactAtomicPPGF}
Construct the exact atomic pseudo-particle Green's function.
Parameters
β
: Inverse temperature.ed
: Exact diagonalization structure describing the atomic problem.
atomic_ppgf(
grid::Keldysh.FullTimeGrid,
ed::KeldyshED.EDCore
) -> Vector{Keldysh.GenericTimeGF{ComplexF64, false, Keldysh.FullTimeGrid}}
Compute atomic pseudo-particle Green's function on a full contour time grid
for a time-independent exact diagonalization problem ed
.
As the resulting PPGF $P(z, z')$ is defined up to a multiplier $e^{-i\lambda (z-z')}$, we choose the energy shift $\lambda$ to fulfil the normalization property $\mathrm{Tr}[i P(-i\beta, 0)] = 1$.
atomic_ppgf(
grid::Keldysh.ImaginaryTimeGrid,
ed::KeldyshED.EDCore
) -> Vector{Keldysh.ImaginaryTimeGF{ComplexF64, false}}
Compute atomic pseudo-particle Green's function on an imaginary time grid
for a time-independent exact diagonalization problem ed
.
As the resulting PPGF $P(\tau)$ is defined up to a multiplier $e^{-\lambda\tau}$, we choose the energy shift $\lambda$ to fulfil the normalization property $\mathrm{Tr}[i P(\beta)] = 1$.
QInchworm.ppgf.partition_function
— Functionpartition_function(
P_0::Vector{QInchworm.exact_atomic_ppgf.ExactAtomicPPGF}
) -> ComplexF64
Extract the partition function $Z = \mathrm{Tr}[i P_0(-i\beta, 0)]$ from a un-normalized pseudo-particle Green's function P_0
.
partition_function(
P::Vector{<:Keldysh.AbstractTimeGF}
) -> ComplexF64
Extract the partition function $Z = \mathrm{Tr}[i P(-i\beta, 0)]$ from a un-normalized pseudo-particle Green's function P
.
QInchworm.ppgf.density_matrix
— Functiondensity_matrix(
P::Vector{QInchworm.exact_atomic_ppgf.ExactAtomicPPGF}
) -> Vector{Matrix{ComplexF64}}
Extract the equilibrium density matrix $\rho = i P(-i\beta, 0)$ from a normalized pseudo-particle Green's function P
. The density matrix is block-diagonal and is returned as a vector of blocks.
density_matrix(
P::Vector{<:Keldysh.AbstractTimeGF}
) -> Vector{Matrix{ComplexF64}}
Extract the equilibrium density matrix $\rho = i P(-i\beta, 0)$ from a normalized pseudo-particle Green's function P
. The density matrix is block-diagonal and is returned as a vector of blocks.
QInchworm.ppgf.normalize!
— Functionnormalize!(
P::Vector{<:Keldysh.AbstractTimeGF},
β::Float64
) -> ComplexF64
Normalize a pseudo-particle Green's function P
by multiplying it by $e^{-i\lambda (z-z')}$ with $\lambda$ chosen such that $\mathrm{Tr}[i P(-i\beta, 0)] = 1$.
Returns
The energy shift $\lambda$.
normalize!(
P::Vector{<:Keldysh.AbstractTimeGF},
τ::Keldysh.TimeGridPoint
)
Normalize a pseudo-particle Green's function P
by multiplying it by $e^{-i\lambda (z-z')}$ with $\lambda$ chosen such that $\mathrm{max}[i P(-i\tau, 0)] = 1$.
normalize!(P_s::Keldysh.AbstractTimeGF, λ)
Multiply a given diagonal block of a pseudo-particle Green's function P_s
by $e^{-i\lambda (z-z')}$.
normalize!(
P_s::QInchworm.spline_gf.IncSplineImaginaryTimeGF{ComplexF64, false},
λ
)
Multiply a given diagonal block of a pseudo-particle Green's function P_s
by $e^{-i\lambda (z-z')}$. This method is defined for the spline-interpolated imaginary-time propagators.
QInchworm.spline_gf
QInchworm.spline_gf
— ModuleQInchworm.spline_gf.SplineInterpolatedGF
— Typestruct SplineInterpolatedGF{GFType, T, scalar} <: Keldysh.AbstractTimeGF{T, scalar}
Wrapper around a Green's function object that allows for fast cubic spline interpolation on the time grid.
The wrapper supports square bracket access to the wrapped object, direct access to the grid
property, eltype()
, Keldysh.norbitals()
and Keldysh.TimeDomain()
. Evaluation at an arbitrary contour time point (via operator ()
) is carried out by a stored set of pre-computed B-spline interpolants.
Fields
GF::Any
: Wrapped Green's functioninterpolants::Any
: B-spline interpolants, one object per matrix element of G
QInchworm.spline_gf.SplineInterpolatedGF
— MethodSplineInterpolatedGF(
GF::Keldysh.AbstractTimeGF{T<:Number, scalar};
τ_max
) -> QInchworm.spline_gf.SplineInterpolatedGF{_A, T} where {_A, T<:Number}
Make a SplineInterpolatedGF
wrapper around GF
and compute interpolants of its data from the start of the grid up to τ_max
. By default, the entire data array is used.
QInchworm.keldysh_dlr
QInchworm.keldysh_dlr
— ModuleExtension of Keldysh.jl defining imaginary time Green's functions represented using the Discrete Lehmann Representation as implemented in Lehmann.jl.
Exports
QInchworm.keldysh_dlr.DLRImaginaryTimeGrid
— Typestruct DLRImaginaryTimeGrid <: Keldysh.AbstractTimeGrid
Wrapper around Lehmann.jl describing a Discrete Lehmann Representation imaginary time grid conforming to the interface of Keldysh.jl TimeGrids.
Fields
contour::Keldysh.ImaginaryContour
points::Vector{Keldysh.TimeGridPoint}
branch_bounds::Tuple{Pair{Keldysh.TimeGridPoint, Keldysh.TimeGridPoint}}
ntau::Int64
dlr::Lehmann.DLRGrid
QInchworm.keldysh_dlr.DLRImaginaryTimeGF
— Typestruct DLRImaginaryTimeGF{T, scalar} <: Keldysh.AbstractTimeGF{T, scalar}
Wrapper around Lehmann.jl describing a Discrete Lehmann Representation imaginary time Green's function conforming to the interface of Keldysh.jl AbstractTimeGF.
Fields
grid::QInchworm.keldysh_dlr.DLRImaginaryTimeGrid
mat::Keldysh.PeriodicStorage
ξ::Keldysh.GFSignEnum
QInchworm.utility
QInchworm.utility
— ModuleQInchworm.utility.ph_conj
— Functionph_conj(
G_int::QInchworm.spline_gf.SplineInterpolatedGF{GFType<:Keldysh.ImaginaryTimeGF{T<:Number, true}, T<:Number, true}
)
Given a spline-interpolated scalar-valued imaginary time Green's function $g(\tau)$, return its particle-hole conjugate $g(\beta-\tau)$.
ph_conj(g::Keldysh.ImaginaryTimeGF{T, true}) -> Any
Given a scalar-valued imaginary time Green's function $g(\tau)$, return its particle-hole conjugate $g(\beta-\tau)$.
ph_conj(
g::QInchworm.keldysh_dlr.DLRImaginaryTimeGF{T, true}
) -> QInchworm.keldysh_dlr.DLRImaginaryTimeGF{T, true} where T
Given a scalar-valued imaginary time Green's function $g(\tau)$, return its particle-hole conjugate $g(\beta-\tau)$.