jED
Documentation for jED.
jED.AIM
jED.AIMParams
jED.Basis
jED.Basis
jED.Eigenspace
jED.Eigenspace
jED.Fockstate
jED.Model
jED.Operator
jED.Overlap
jED.CDag_sign
jED.C_sign
jED.GLoc
jED.GWeiss
jED.GWeiss
jED.GWeiss_from_Imp
jED.N_do
jED.N_el
jED.N_up
jED.S
jED._H_CDagC
jED._H_nn
jED._block_slice
jED._find_cdag_overlap_blocks
jED._generate_blocks!
jED._overlap_cdagger_ev!
jED._overlap_list
jED.ann
jED.ann_op
jED.calc_E
jED.calc_GF_1
jED.calc_Hamiltonian
jED.calc_Z
jED.create
jED.create_op
jED.operator_ni
jED.overlap
jED.overlap_2
jED.overlap_cdagger_c
jED.overlap_ni_nj
jED.Δ_AIM
jED.Δ_from_GWeiss
jED.Σ_from_GImp
jED.AIM
— TypeAIM <: Model
AIM(ϵₖ::Vector{Float64}, Vₖ::Vector{T}, μ::Float64, U::Float64)
Type for the Anderson impurity model. A model can be used to construct a Eigenspace
given a set of Basis
.
Can be created from
- bath onsite energies $\epsilon_k$,
- bath-impurity hopping $V_k$
- impurity chemical potential $\mu$
- Coulomb interaction strength on the impurity $U$
Example
julia> ϵₖ = [0.5, -5.0]
julia> Vₖ = [1.0, 1.0]
julia> U = 1.0
julia> μ = 0.5
julia> model = AIM(ϵₖ, Vₖ, μ, U)
AIM{3, Float64}([
-0.5 1.0 1.0; 1.0 0.5 0.0; 1.0 0.0 -5.0], [1.0 0.0 0.0; 0.0 0.0 0.0; 0.0 0.0 0.0], [0.0 0.0 0.0; 0.0 0.0 0.0; 0.0 0.0 0.0])
jED.AIMParams
— TypeAIMParams
Parameters for energy levels ($\epsilon_k$) and hybridization amplitudes ($V_k$) of the Anderson impurity model.
jED.Basis
— TypeBasis
A state is represented as a BitVetor of length 64. The first 32 bit (ordering follows Julia Vector indexing) correspond to spin up for each site (index 1 corresponds to the first site, 2 to the second, etc.) and bits 33 to 64 to spin down, for site 1 to 32. Each bit is set to 1 if it is occupied by an electron of the given spin. In principle, the number of electron flavors and orbitals is free, however, for now it is fixed to 2 and 1.
A list of states should be generated using States(NSites)
with NSites being the number of sites in the model.
Fields
NFlavors
: Integer, number of flavors (TODO: hardcoded to2
for now. later: orbitals*2?!)NSites
: Integer, number of sitesstates
: Vector{State}, list of Fock states. The are sorted according to good quantum numbers.blocklist
contains start and size of blocks.blocklist
: Vector{NTuple{4,Int}}, Vector of 4-tuples. Each entry encodes a block in the following form: 1. element is the start index of the block, 2. element is the column size of the block, 3. element is the electron number (seeN_el
), 4. element is the spin (seeS
)cdag_ov
: Vector{Tuple{UnitRange{Int},UnitRange{Int}}}, Slices of states, that contribute to an overlap with one creation operator. Needed for performance reasons in computation of Lehmann representation (see also (@ref findcdagoverlapblocks)[_find_cdag_overlap_blocks
])).
jED.Basis
— MethodBasis(NSites::Int; NFlavors::Int = 2)
Contructs a Fock basis for NSites
with NFlavors
.
jED.Eigenspace
— TypeEigenspace
Containes Eigenvalues and Eigenvectors grouped into blocks with are each ordered by magnitude of Eigenvalues. The blocks are part of the Basis
.
Fields
evals
: Eigenvaluesevecs
: EigenvectorsE0
: smallest Eigenvalue
jED.Eigenspace
— MethodEigenspace(model::Model, basis::Basis; verbose::Bool = true, FPT::Type{FPTi} = eltype(model.tMatrix))
Constructs the Eigenspace for Model
over given Basis
by diagonalizing the Hamiltonian (see also calc_Hamiltonian
) for each block.
jED.Fockstate
— TypeFockstate
Internal representation of a Fockstate.
jED.Model
— TypejED.Operator
— TypeOperator
Operator for Fockstates
. Holds function for transformation of Fock state and change in quantum numbers.
Fields
f
: Function, operation on Fock stateN_inc
: Integer, change in electron numberS_inc
: Integer, change in spin number
Example
julia> o1 = Operator(x->create())
jED.Overlap
— TypeOverlap
Contains information about overlap of Fock states under given Operator
.
Fields
op
: Function, operation on Fock stateov_blocks
: Vector{Int}, index is starting block, entry targetov_list
: Vector{Int}, index is starting block, entry index in eigenvector (seeEigenspace
)
jED.CDag_sign
— MethodCDag_sign(state,i)
Sign when annihilating an electron at position i
in state
. Returns: -1/1/0 for Uneven/Even permutations and 0 when the state is not occupied.
julia> CDag_sign(Fockstate{6}(Bool[1,0,1,1,0,0]),5)
-1
julia> CDag_sign(Fockstate{6}(Bool[1,0,1,1,0,0]),3)
0
julia> CDag_sign(Fockstate{7}(Bool[1,0,1,0,1,0,0]),4)
1
jED.C_sign
— MethodC_sign(state,i)
Sign when creating an electron at position i
in state
. Returns: -1/1/0 for Uneven/Even permutations and 0 when the state is not occupied.
julia> C_sign(Fockstate{6}(Bool[1,0,1,1,0,0]),3)
-1
julia> C_sign(Fockstate{6}(Bool[1,0,1,1,0,0]),2)
0
julia> C_sign(Fockstate{6}(Bool[1,0,1,1,0,0]),4)
jED.GLoc
— MethodGLoc(ΣImp::MatsubaraF, μ::Float64, νnGrid::FermionicMatsubaraGrid, kG::KGrid)
Compute local Green's function $\int dk [i\nu_n + \mu + \epsilon_k - \Sigma_\text{Imp}(i\nu_n)]^{-1}$. TODO: simplify -conj!!!
jED.GWeiss
— MethodGWeiss(νnGrid::FermionicMatsubaraGrid, p::AIMParams)
GWeiss!(target::Vector, νnGrid::Vector, μ::Float64, p::AIMParams)
Computes Weiss Green's frunction from Anderson Parameters
.
jED.GWeiss
— MethodGWeiss(Δ::OffsetVector{ComplexF64, Vector{ComplexF64}}, μ::Float64, νnGrid::FermionicMatsubaraGrid)
Computes Weiss Green's frunction from hybridization function
.
jED.GWeiss_from_Imp
— MethodGWeiss_from_Imp(GLoc::MatsubaraF, ΣImp::MatsubaraF)
Compute Updates Weiss Green's function from impurity self-energy via $[G_\text{loc} + \Sigma_\text{Imp}]^{-1}$ (see GLoc
and Σ_from_GImp
).
jED.N_do
— MethodN_do(s::Fockstate{NSites})
Number of down electrons, $N_\downarrow$ in state s
.
jED.N_el
— MethodN_el(s::AbstractVector)::Int
Total electron number of state.
jED.N_up
— MethodN_up(s::Fockstate{NSites})
Number of up electrons, $N_\uparrow$ in state s
.
jED.S
— MethodS(s::AbstractVector)::Int
Total spin of state. #TODO: only implemented for flavor=2
jED._H_CDagC
— Method_H_CDag_C(istate, jstate, tMatrix)
Returns the hopping contribution for states $\sum_{i,j} \langle i | T | j \rangle$, with T being the hopping matrix tmatrix
and the states i and j given by istate
and jstate
.
jED._H_nn
— Method_H_nn(istate, jstate, tMatrix)
Returns the hopping contribution for states $\sum_{i,j} \langle i | T | j \rangle$, with T being the hopping matrix tmatrix
and the states i and j given by istate
and jstate
.
jED._block_slice
— Method_block_slice(bi::Blockinfo)
Slice of continuous vector for given block bi
.
jED._find_cdag_overlap_blocks
— Method_find_cdag_overlap_blocks(blocklist::Vector{Blockinfo}; S_inc=-1, N_inc = 1)
Find block index with N_el
increased by 1
and S
either in increased or decreased (depending on S_inc
), e.g. block index of state $\langle j |$ with non-zero overlap of $c^\dagger_\uparrow |i \rangle$.
Returns: Vector{Int}, equal length to blocklist
. Each entry with index i
contains the index j
of the block with one more electron and spin (i.e. the block for which $\langle j| c^\dagger | i \rangle$ does NOT vanish). This is stored here for performance reasons, since these overlaps are used often in the Lehman representation. TODO: do not hardcode spin up GF!!! (this forces S+1 instead of S+-1 for now)
jED._generate_blocks!
— Method_generate_blocks!(states::Vector{Fockstate})
Sort state list and generate list of blocks.
jED._overlap_cdagger_ev!
— Method_overlap_cdagger_ev!(tmp::Vector{Float64}, vec1::Vector{Float64}, vec2::Vector{Float64}, perm::Vector{Int})::Float64
Computes $\langle \text{EV}_j | c^\dagger_k | \text{EV}_i \rangle$ where $\c^\dagger$ has been computed in the Fockbasis and is given as a permutation of indices for $\langle \text{EV}_j |$.
This is an internal function used to compute the overlaps
.
jED._overlap_list
— Methodoverlaplist(basis::Basis[, i::Int, j::Int], op::Function)
Calculates the overlap between block i
and j
under op
. Builds full overlap indices for all block if i
and j
are not provided.
Example
julia> basis = jED.Basis(3)
julia> op(b) = jED.create(b, 4) # creation operator for spin down at site 1 of 3
julia> jED._overlap_list(basis, 2, 4, op)
[0, 1, 2]
jED.ann
— Methodann(state::Fockstate{Length}, i::Int)::Union{Nothing,Fockstate} where Length
Returns either a new Fockstate
with an electron i
annihilated (see layout of Fockstate, how up and down is encoded), or nothing
, if i
is not occupied. See also create
. See also ann_op
for the Operator
version.
jED.ann_op
— Methodanne_op(b::Basis, i::Int)
Operator
for annihilation of state i
(for basis length N, i > N annihilates spin down, otherwise spin up).
jED.calc_E
— Methodcalc_E(eigenspace, β)
Calculates the total energy from eigenspace
and the inverse temperature β
.
jED.calc_GF_1
— MethodcalcGF1(basis::Basis, es::Eigenspace, νnGrid::AbstractVector{ComplexF64}, β::Float64; ϵ_cut::Float64=1e-16, overlap=nothing)
Computes $|\langle i | c^\dagger | j \rangle|^2 \frac{e^{-\beta E_i} + e^{-\beta E_j}}{Z (E_j - E_i + freq)}$.
Arguments
basis
: Basis, obtained withBasis
.es
: Eigenspace, obtained withEigenspace
νnGrid
: AbstractVector{ComplexF64}, list of Matsubara Frequenciesβ
: Float64, inverse temperature.ϵ_cut
: Float64, cutoff for $e^{-\beta E_n}$ terms in Lehrmann representation (all contributions below this threshold are disregarded)overlap
: Overlap, precalculated overlap between blocks of basis. Obtained withOverlap
jED.calc_Hamiltonian
— Methodcalc_Hamiltonian(model::Model, basis::Basis)
Calculates the Hamiltonian for a given
jED.calc_Z
— Methodcalc_Z(eigenspace, β)
Calculates the partition function from eigenspace
and the inverse temperature β
.
jED.create
— Methodcreate(state::Fockstate{Length}, i::Int)::Union{Nothing,Fockstate} where Length
Returns either a new Fockstate
with an electron i
created (see layout of Fockstate, how up and down is encoded), or nothing
, if i
is already occupied. See also ann
. See also create_op
for the Operator
version.
jED.create_op
— Methodcreate_op(b::Basis, i::Int)
Operator
for creation of state i
(for basis length N, i > N creates spin down, otherwise spin up).
jED.operator_ni
— Methodoperator_ni(state::Fockstate, i::Int)
Calculate n_i |ket⟩, i.e. i
is the index for the number operator
Returns: (eigenval,newState)::Tuple{Int, Fockstate}
with eigenval
1
(state i
occupied) or 0
(state i
unoccupied) and newState
. In case of this operator, newState === state
. The density overlap of 2 states can be calculated efficiently using overlap_ni_nj
.
Example
julia> jED.operator_ni(jED.SVector{8}(Bool[1,1,1,1,0,1,0,1]), 2);
(1, Bool[1, 1, 1, 1, 0, 1, 0, 1])
jED.overlap
— Methodoverlap(bra::Fockstate, ket::Fockstate)
Calculates ⟨bra
|ket
⟩. This is useful when operators have been applied to one state and no efficient direct implementation for an overlap involving this operator is available.
Returns: True/False
jED.overlap_2
— Methodoverlap_2(bra::Vector, ket::Vector)
Computes |⟨bra
|ket
⟩|² for Float64
or ComplexF64
type vectors.
jED.overlap_cdagger_c
— Methodoverlap_cdagger_c(bra::Fockstate, i::Int, ket::Fockstate, j::Int)::Int
Calculates $\langle$ bra
$| c^\dagger_i c_j |$ ket
$\rangle$, i.e. i
is the index for the creation operator and j
the index for the annihilation operator. Internally, we check, that both states have exactly two or no occupation difference. In both cases the product of C_sign
and CDag_sign
is returned, otherwise 0. Note, that this allows spin flips!
Returns: -1/0/1
Example
julia> t1 = jED.SVector{8}(Bool[1,1,1,1,0,1,0,1]);
julia> t2 = jED.SVector{8}(Bool[1,1,0,1,1,1,0,1]);
julia> jED.overlap_cdagger_c(t1,t2,3,3)
0
julia> jED.overlap_cdagger_c(t1,t1,2,2)
1
jED.overlap_ni_nj
— Methodoverlap_ni_nj(bra::Fockstate, ket::Fockstate, i::Int, j::Int)
Calculate ⟨bra| n^†i nj |ket⟩
Returns: True/False (converts to 1/0)
jED.Δ_AIM
— MethodΔ_AIM(νnGrid::FermionicMatsubaraGrid, p::AIMParams)
Δ_AIM(νnGrid::FermionicMatsubaraGrid, p::Vector{Float64})
Computes hybridization function $\sum_p \frac{V_p^2}{i\nu_n - \epsilon_p}$ from Anderson Impurity Model Parameters
with $p$ bath sites.
jED.Δ_from_GWeiss
— MethodΔ_from_GWeiss(GWeiss::MatsubaraF, μ::Float64, νnGrid::FermionicMatsubaraGrid)
Computes hybridization function from Weiss Green's function via $\Delta^{\nu} = i\nu_n + \mu - \left(\mathcal{G}^{\nu}_0\right^{-1}$.
jED.Σ_from_GImp
— MethodΣ_from_GImp(GWeiss::OffsetVector{ComplexF64, Vector{ComplexF64}}, GImp::OffsetVector{ComplexF64, Vector{ComplexF64}})
Computes self-energy from impurity Green's function (as obtained from a given impurity solver) and the Weiss Greens Function
.