RandomPhaseApproximation

Documentation for RandomPhaseApproximation.

Introduction

Standard random phase approximation (particle-hole channel) for quantum lattice systems based on the QuantumLattices and TightBindingApproximation packages.

Installation

In Julia v1.8+, please type ] in the REPL to use the package mode, then type this command:

pkg> add https://github.com/Quantum-Many-Body/RandomPhaseApproximation.jl

Getting Started

Examples of random phase approximation for quantum lattice system

Manuals

RandomPhaseApproximation.EigenRPAMethod
EigenRPA(path::Union{ReciprocalPath, ReciprocalZone}, bz::ReciprocalZone; onlyvalue::Bool=true,  options...)

Construct a EigenRPA type. Attribute options contains (gauge=:icoordinate, exchange=false, η=1e-8, temperature=1e-12, μ=0.0, bnd=nothing)

source
RandomPhaseApproximation.PHVertexRepresentationType
PHVertexRepresentation{H<:RepresentationGenerator, Vq, Vk, T} <: MatrixRepresentation

Matrix representation of the particle-hole channel of two-body interaction terms. When the k₁ and k₂ are nothing, the exchange terms is ommitted. $1/N\sum_{kk'q\alphaeta m n}[V^{ph}_{\alpha\beta,mn}(q)-V^{ph}_{\alpha m,\beta n}(k'-k)]c^\dagger_{k-q\alpha}c_{k\beta}c^\dagger_{k' n}c_{k'-q m}$ c^†ᵢ = 1/√N∑ₖc†ₖ exp(-ik*rᵢ)

source
RandomPhaseApproximation.PHVertexRepresentationMethod
PHVertexRepresentation{H}(table, gauge::Symbol=:icoordinate)  where {H<:RepresentationGenerator}
PHVertexRepresentation{H}(q, table, gauge::Symbol=:icoordinate) where {H<:RepresentationGenerator}

Get the matrix representation of particle-hole channel.

source
RandomPhaseApproximation.ParticleHoleSusceptibilityType
ParticleHoleSusceptibility{P<:Union{ReciprocalPath,ReciprocalZone}, RZ<:ReciprocalZone, E<:AbstractVector, S<:Operators} <: Action

Calculate the particle-hole susceptibility within random phase approximation. Attribute options contains (η=0.01, gauge =:icoordinate, temperature=1e-12, μ=0.0, findk = false)

source
RandomPhaseApproximation.ParticleHoleSusceptibilityMethod
ParticleHoleSusceptibility(path::Union{ReciprocalPath,ReciprocalZone}, bz::ReciprocalZone, energies::AbstractVector, operators::Tuple{AbstractVector{<:Operators}, AbstractVector{<:Operators}}; options...)

Construct a ParticleHoleSusceptibility type.

source
RandomPhaseApproximation.RPAMethod
RPA(
    tba::AbstractTBA{K, <:OperatorGenerator}, 
    uterms::Tuple{Vararg{Term}}
) where {K<:TBAKind}
RPA(
    lattice::AbstractLattice, 
    hilbert::Hilbert, 
    terms::Tuple{Vararg{Term}}, 
    uterms::Tuple{Vararg{Term}}; 
    neighbors::Union{Nothing, Int, Neighbors}=nothing, 
    boundary::Boundary=plain
)   
RPA(
    tba::AbstractTBA{K, <:AnalyticalExpression}, 
    hilbert::Hilbert, 
    table::Table, 
    uterms::Tuple{Vararg{Term}}; 
    neighbors::Union{Nothing, Int, Neighbors}=nothing, 
    boundary::Boundary=plain
) where {K<:TBAKind}

Construct a RPA type.

source
QuantumLattices.add!Method
add!(dest::AbstractMatrix, mr::PHVertexRepresentation{<:RepresentationGenerator}, m::Operator; kwargs...)

Get the matrix representation of an operator and add it to destination.

source
RandomPhaseApproximation._chikq0Method
_chikq0(tba::AbstractTBA, k::AbstractVector, q::AbstractVector, omega::Float64; 
eta::Float64=0.01, tem::Float64=1e-12, mu::Float64=0.0, kwargs...)

Return chi0(k,q){ij,mn}= chi0(k,q){-+,+-}, chi0(k,q){12,34}== <c^\dagger{k,2}c{k-q ,1}c^\dagger{k-q,3}c_{k,4}>

source
RandomPhaseApproximation.chikqmMethod
chikqm(tba::AbstractTBA, bz::ReciprocalZone, path::Union{ReciprocalZone, ReciprocalPath}, vph::Array{T,5}; tem::Float64=1e-12, mu::Float64=0.0, onlyvalue::Bool=true, η::Float64=1e-6, bnd::Union{UnitRange{Int}, StepRange{Int,Int}, Vector{Int}, Nothing}=nothing, kwargs...) where T<:Number -> Tuple{ Array{Array{ComplexF64, 1}, 1}, Array{Matrix{ComplexF64}, 1}, Array{Matrix{ComplexF64}, 1}}

Get the eigenvalues, eigenvectors, and unitary transformation (from orbital to band) of particle-hole susceptibilities χ_{αβ}(ω,k1,k2,q). Now only the zero-temperature case is supported.
vph store the particle-hole vertex,e.g. vph[ndim,ndim,nk,nk,nq] where sqrt(nidm) is the number of degrees of freedom in the unit cell,nk=length(bz), nq=length(path).
tem is the temperature;
mu is the chemical potential.
onlyvalue : only need the eigenvalues ( isval=true,zero temperature ), isval=false denotes that the cholesky method is used.
η:small number to advoid the semipositive Hamiltonian,i.e. Hamiltonian+diagm([η,η,...])
bnd:select the bands to calculate the χ₀
kwargs store the keys which are transfered to matrix function.

source
RandomPhaseApproximation.chiqMethod
chiq(eigenvc::Array{ComplexF64, 3}, eigenval::Array{Float64, 2}, bz::ReciprocalZone, path::Union{ReciprocalPath, ReciprocalZone}, vph::Array{T, 3}, omegam::Vector{Float64}; eta::Float64=0.01, tem::Float64=1e-12, mu::Float64=0.0, scflag=false) where T<:Number -> Tuple{ Array{ComplexF64, 4}, Array{ComplexF64, 4} }
source
RandomPhaseApproximation.chiqMethod
chiq(vph::Array{T, 3}, chi0::Array{ComplexF64, 4}) where T<:Number -> Array{ComplexF64, 4}

Get the susceptibility χ_{αβ}(ω,q).

source
RandomPhaseApproximation.chiqMethod
chiq(
    tba::AbstractTBA, 
    bz::AbstractVector{<:AbstractVector}, 
    path::Union{ReciprocalPath, ReciprocalZone}, 
    vph::Array{T,3}, omegam::AbstractVector; 
    eta::Float64=0.01, 
    tem::Float64=1e-12, 
    mu::Float64=0.0, 
    scflag=false, 
    kwargs...
) where T<:Number -> Tuple{ Array{ComplexF64, 4}, Array{ComplexF64, 4} }

Get the particle-hole susceptibilities χ⁰(ω,q) and χ(ω,q). The spectrum function is satisfied by A(ω,q) = Im[χ(ω+i*0⁺,q)].

Arguments

  • vph is the bare particle-hole vertex (vph[ndim,ndim,nq])
  • omegam store the energy points
  • eta is the magnitude of broaden
  • tem is the temperature
  • mu is the chemical potential
  • scflag == false (default, no superconductivity), or true ( BdG model)
  • kwargs is transfered to matrix(tba; kwargs...) function
source
RandomPhaseApproximation.chiq0Method
chiq0(tba::AbstractTBA, bz::ReciprocalZone, path::Union{ReciprocalPath, ReciprocalZone}, omegam::Vector{Float64}; eta::Float64=0.01, tem::Float64=1e-12, mu::Float64=0.0, scflag=false, kwargs...) -> Array{Float64, 4}

Get the particle-hole susceptibilities χ⁰(ω,q).
omegam store the energy points;
eta is the magnitude of broaden;
tem is the temperature;
mu is the chemical potential.
'kwargs' is transfered to matrix(tba;kwargs...) function
scflag == false(default) => particle-hole channel, ==true => Nambu space.

source
RandomPhaseApproximation.correlationMethod
correlation(χ::Array{<:Number, 4}, path::Union{ReciprocalPath, ReciprocalZone}, operators::Tuple{AbstractVector{S}, AbstractVector{S}}, table; gauge=:rcoordinate) where {S <: Operators} -> Matrix

Return physical particle-hole susceptibility.

source
RandomPhaseApproximation.findkMethod
findk(kq::AbstractVector, bz::ReciprocalZone) -> Int

Find the index of k point in the reduced Brillouin Zone bz, i.e. bz[result] ≈ kq

source
RandomPhaseApproximation.projchiMethod
projchi(val::Array{Array{T, 1}, 1}, vec::Array{Array{T1, 2}, 1}, f2c::Array{ Array{ComplexF64, 2}, 1}, omegam::Vector{Float64}, eta::Float64=1e-2) where {T<:Number,T1<:Number} -> Array{ComplexF64,4}

Get the $\chi_{ij,nm}(\omega,q)$. The eigenvalues, eigenvectors, and unitary matrix are obtained by methodchikqm`.

source
RandomPhaseApproximation.projchiimMethod
projchiim(val::Array{Array{T, 1}, 1}, vec::Array{Array{T1, 2}, 1}, f2c::Array{ Array{ComplexF64, 2}, 1},omegam::Vector{Float64}, eta::Float64=1e-2) where {T<:Number, T1<:Number} -> Array{ComplexF64, 4}

Get the imaginary party of $\chi_{ij,nm}(\omega,q)$. The eigenvalues,eigenvectors, and unitary matrix are obtained by methodchikqm`.

source
RandomPhaseApproximation.selectpathMethod
selectpath(path::AbstractVector{<:Tuple{<:AbstractVector, <:AbstractVector}}, bz::ReciprocalZone;ends::Union{<:AbstractVector{Bool},Nothing}=nothing, atol::Real=atol, rtol::Real=rtol) -> Tuple(Vector{Vector{Float64}}, Vector{Int})  -> -> Tuple(ReciprocalPath, Vector{Int})

Select a path from the reciprocal zone. Return ReciprocalPath and positions of points in the reciprocal zone.

source
RandomPhaseApproximation.selectpathMethod
selectpath(stsp::Tuple{<:AbstractVector, <:AbstractVector}, rz::ReciprocalZone; ends::Tuple{Bool, Bool}=(true, false), atol::Real=atol, rtol::Real=rtol) -> Tuple(Vector{Vector{Float64}}, Vector{Int})

Select a path from the reciprocal zone.

source
RandomPhaseApproximation.vertex_phFunction
vertex_ph(rpa::RPA, path::AbstractVector{<:AbstractVector}, gauge=:icoordinate) -> Array{ComplexF64, 3}

Return particle-hole vertex induced by the direct channel of interaction ( except the Hubbard interaction which include direct and exchange channel).

source
RecipesBase.apply_recipeFunction
@recipe plot(pack::Tuple{Algorithm{<:RPA}, Assignment{<:ParticleHoleSusceptibility}}, mode::Symbol=:χ)

Define the recipe for the visualization of particle-hole susceptibilities.

source
RecipesBase.apply_recipeFunction
@recipe plot(pack::Tuple{Algorithm{<:RPA}, Assignment{<:EigenRPA}})

Define the recipe for the visualization of particle-hole susceptibilities.

source
RecipesBase.apply_recipeFunction
@recipe plot(rz::ReciprocalZone, path::Union{ReciprocalPath,Nothing}=nothing)

Define the recipe for the visualization of a reciprocal zone and a reciprocal path

source