SUNSpinWaveTheory
Documentation for SUNSpinWaveTheory.
Introduction
SU(N) spin wave theory for free quantum lattice systems based on the QuantumLattices pack and the TightBindingApproximation pack.
Installation
Getting Started
Examples of SU(N) spin-wave theory for quantum lattice system
Manuals
QuantumLattices.DegreesOfFreedom.Hilbert
— MethodHilbert(hilbert::Hilbert{<:Fock{:b}}, magneticstructure::MagneticStructure)
Get the corresponding Hilbert space of the original one after the Holstein-Primakoff transformation.
QuantumLattices.DegreesOfFreedom.Metric
— MethodMetric(::SUNMagnonic, hilbert::Hilbert{<:Fock{:b}}) -> OperatorUnitToTuple
Get the index-to-tuple metric for a quantum spin system after the Holstein-Primakoff transformation.
SUNSpinWaveTheory.HPTransformation
— TypeHPTransformation{S<:Operators, U<:CompositeIndex, M<:MagneticStructure} <: UnitSubstitution{U, S}
Holstein-Primakoff transformation.
SUNSpinWaveTheory.InelasticNeutron
— TypeInelasticNeutron <: SpectraKind{:INS}
Inelastic Neutron Scattering Spectra of quantum lattice system.
SUNSpinWaveTheory.MagneticStructure
— TypeMagneticStructure{L<:AbstractLattice, P<:Int, D<:Number, T<:Complex}
The magnetic structure of an ordered quantum lattice system. The cell
is a magnetic cell.
SUNSpinWaveTheory.MagneticStructure
— MethodMagneticStructure(cell::AbstractLattice, moments::Dict{<:Int, <:AbstractVector})
(ms::MagneticStructure)(orders::Dict{<:Int, <:AbstractMatrix}) -> Dict{<:Int, <:Number}
1).Construct the magnetic structure on a given lattice with the given moments. The moment is given by the angles (θ₁,θ₂,...,ψ₁,ψ₂,...). 2).Get order parameters for every site. e.g. |gs⟩ = ms.rotations[i][1,:] where gs = classical ground state; ⟨σᵢˣ⟩ = ⟨gs|σᵢˣ|gs⟩.
SUNSpinWaveTheory.Multipole
— TypeMultipole <: SpectraKind{:Multipole}
Multipole Scattering Spectra of quantum lattice system. (∑{α}s{α}*s_{α})
SUNSpinWaveTheory.SUNLSWT
— TypeSUNLSWT{K<:TBAKind{:BdG}, L<:AbstractLattice, Hₛ<:OperatorGenerator, HP<:HPTransformation, Ω<:Image, H<:Image} <: AbstractTBA{K, H, AbstractMatrix}
SU(N) Linear spin wave theory for magnetically ordered quantum lattice systems.
SUNSpinWaveTheory.SUNLSWT
— MethodSUNLSWT(
lattice::AbstractLattice,
hilbert::Hilbert,
terms::Tuple{Vararg{Term}},
magneticstructure::MagneticStructure;
neighbors::Union{Nothing, Int, Neighbors}=nothing,
boundary::Boundary=plain
)
Construct a SUNLSWT. lattice
is the original lattice.
SUNSpinWaveTheory.SUNMagnonic
— TypeSUNMagnonic <: TBAKind{:BdG}
Magnonic quantum lattice system.
SUNSpinWaveTheory.Spectra
— TypeSpectra{K<:SpectraKind, P<:ReciprocalSpace, E<:AbstractVector, S<:Operators, O} <: Action
Spectra of 'magnetically' ordered quantum lattice systems by SU(N) linear spin wave theory.
SUNSpinWaveTheory.SpectraKind
— TypeSpectraKind{K}
The kind of a spectrum calculation.
Base.nameof
— Methodnameof(alg::Algorithm, assign::Assignment, ecut::Float64) -> String
Get the name of the combination of an algorithm, an assignment and ecut.
QuantumLattices.add!
— Methodadd!(dest::Matrix,
mr::TBAMatrixRepresentation{SUNMagnonic},
m::Operator{<:Number, <:ID{CompositeIndex{<:Index{Int, <:FID{:b}}}, 2}};
atol=atol/5,
kwargs...
) -> typeof(dest)
Get the matrix representation of an operator and add it to destination.
RecipesBase.apply_recipe
— Function@recipe plot(pack::Tuple{Algorithm{<:SUNLSWT}, Assignment{<:Spectra}}, Ecut::Float64, dE::Float64=1e-3)
Define the recipe for the visualization of an assignment with Ecut and data of an algorithm.
RecipesBase.apply_recipe
— Method@recipe plot(pack::Tuple{Algorithm{<:SUNLSWT}, Assignment{<:Spectra}}, parameters::AbstractVector{Float64})
Define the recipe for the visualization of an assignment with parameters (title) of an algorithm.
SUNSpinWaveTheory.SUNTerm
— FunctionSUNTerm(
id::Symbol,
J₁₂₃₄::Array{<:Number, 4},
bondkind::Int=1;
value = 1,
amplitude::Union{Function, Nothing}=nothing,
modulate::Union{Function, Bool}=false
) -> Coulomb
Hamiltionian: $∑(J_{i1,i2;j3,j4}b†_{i1}b_{i2}b†_{j3}b_{j4})$
SUNSpinWaveTheory.SUNTerm
— FunctionSUNTerm(
id::Symbol,
Jᵤᵥ::Array{<:Number, 2},
bondkind::Int=0;
value = 1,
amplitude::Union{Function, Nothing}=nothing,
modulate::Union{Function, Bool}=false
) -> Onsite
SUNSpinWaveTheory.SUNTerm
— FunctionSUNTerm(
id::Symbol,
J₁₂₃₄::Array{<:Number, 2},
dimᵢ::Int,
dimⱼ::Int,
bondkind::Int=1;
value=1,
amplitude::Union{Function, Nothing}=nothing,
modulate::Union{Function, Bool}=false
) -> Term
Return the Term
type. J₁₂₃₄
is the coefficience of Hamiltionian including the exchange term and onsite term, i.e. $∑(J_{i1,j2;i3,j4}b†_{i1}b†_{j2}b_{i3}b_{j4}) + ∑(Bᵢᵤᵥb†ᵢᵤbᵢᵥ )$ dimᵢ
and dimⱼ
are the dimensions of local Hilbert space on sites i and j, respectively.
SUNSpinWaveTheory.fluctuation
— Methodfluctuation(sunlswt::SUNLSWT, bz::BrillouinZone, site::Int, nbonson::Int; imagtol=1e-8) -> Float64
fluctuation(sunlswt::SUNLSWT, bz::BrillouinZone; imagtol=1e-8) -> Dict{<:Int, Float64}
To calculate $⟨b_0^†b_0⟩$ with respect to the quantum ground state |0⟩ at zero temperature.
SUNSpinWaveTheory.kdosOperator
— MethodkdosOperator(u::Matrix{<:Number}, site::Int, lattice::Union{<:AbstractLattice, Int}) -> Operators
Return Operators for kDOS case. u = ⟨α|site,orbital⟩, where |site,orbital⟩ => b_{site,orbital}, |α⟩ is the target state.
SUNSpinWaveTheory.localorder
— Methodlocalorder(sunlswt::SUNLSWT, orders::Dict{Int, <:Matrix{T}}) where T<:Number -> Dict{Int, <:Number}
Return order parameters of each site, <classical gs|o|classical gs>. orders[pid]
is the matrix of physical observable.
SUNSpinWaveTheory.matrixcoef
— Methodmatrixcoef(sunlswt::SUNLSWT) -> Hermitian
Return the coefficience of exchange interactions, i.e. ∑(Jᵢ₁ⱼ₂ⱼ₃ᵢ₄b†ᵢ₁b†ⱼ₂bⱼ₃bᵢ₄) + ∑(Bᵢᵤᵥb†ᵢᵤbᵢᵥ )
SUNSpinWaveTheory.multipoleOperator
— MethodmultipoleOperator(j₁::Union{Int,Rational{Int}}, j₂::Union{Int,Rational{Int}}, j::Union{Int,Rational{Int}}, m::Union{Int,Rational{Int}}) -> Array{Float64}
Construct the matrix of multipole operator, i.e. $M_{jm}^{j₁j₂} = ∑_{m₁m₂} (-1)^{j₂-m₂}C_{m₁,-m₂,m}^{j₁j₂j} b†_{j₁m₁}b_{j₂m₂}$ with $M_{jm}^{†j₁j₂}=M_{j-m}^{j₂j₁}(-1)^{j₂-j₁+m}$. e.g. the basis of dipole matrix is [|1⟩, |0⟩, |-1⟩].
SUNSpinWaveTheory.optimorder
— Methodoptimorder(sunlswt::SUNLSWT; numrand::Int = 0, method = LBFGS(), g_tol = 1e-12, optionskwargs... ) -> Tuple{SUNLSWT,Union{Optim.MultivariateOptimizationResults,Nothing}}
Optimize the ground state of Hamiltionian, i.e. ⟨T|H|T⟩, see DOI: 10.1103/PhysRevB.97.205106 numrand
is the number of the optimization. optionskwargs
is the keywords of Optim.Options(gtol=gtol,optionskwargs...)
References
- Zhao-Yang Dong, Wei Wang, and Jian-Xin Li, SU(N) spin-wave theory: Application to spin-orbital Mott insulators, Physical Review B 97, 205106 (2018)
SUNSpinWaveTheory.optimorder2
— Methodoptimorder2(sunlswt::SUNLSWT, ub::Vector{Float64}, x₀::Vector{Float64}; numrand::Int = 0, rule::Function=x->x, method = LBFGS(), g_tol = 1e-12, optionskwargs... ) -> Tuple{SUNLSWT, Union{Optim.MultivariateOptimizationResults, Nothing}}
Optimize the ground state of Hamiltionian. If the all degrees of freedom is considered, the Function optimorder
is recommended. Function rule
gives the transformation of the selected degrees of freedom to all degrees of freedom.
SUNSpinWaveTheory.optimorder2
— Methodoptimorder2(sunlswt::SUNLSWT; numrand::Int = 0,rule::Function=x->x, method = LBFGS(), g_tol = 1e-12, optionskwargs... ) -> Tuple{SUNLSWT, Union{Optim.MultivariateOptimizationResults, Nothing}}
Optimize the ground state of Hamiltionian, i.e. ⟨T|H|T⟩. If the all degrees of freedom is considered, the Function optimorder
is recommended.
SUNSpinWaveTheory.quantumGSenergy
— MethodquantumGSenergy(sunlswt::SUNLSWT, bz::BrillouinZone; imagtol=1e-8) -> Float64
quantumGSenergy(sunlswt::SUNLSWT, nk::Int; imagtol=1e-8) -> Float64
Calculate the quantum-ground-state energy per magnetic unit cell. $E=E_{calssical}/N_{mc}-1/(2N_{mc})*∑_{kα}ω_{kα} - 1/(4N_{mc})∑_kTr(H(k))$, where $N_{mc}$ is the number of magnetic unit cell
SUNSpinWaveTheory.rotation
— Methodrotation(sitapsi::T, generators::Vector{Matrix{Complex{Int}}}}) where T<:AbstractVector{<:Number} -> Matrix{ComplexF64}
Get the rotation matrix which rotates the |0⟩
state to the target state |T⟩=U(θ₁,θ₂,...,ψ₁,ψ₂,...)b⁰† |0⟩; T†ᵥ = ∑ᵤUᵥᵤ*bᵘ†. The detail is given in doi:10.1103/PhysRevB.97.205106.
SUNSpinWaveTheory.rotation_gen
— Methodrotation_gen(n::Int) -> Vector{Matrix{Complex{Int}}}
Return the generator matrices of SU(n) group.
SUNSpinWaveTheory.spectraEcut
— MethodspectraEcut(ass::Assignment{<:Spectra}, Ecut::Float64, dE::Float64) -> Tuple{Vector{Float64}, Vector{Float64}, Matrix{Float64}}
Construct the spectra with fixed energy. The energy of abs(energy-Ecut) <= dE
is selected. nx
and ny
are the number of x and y segments of ReciprocalZone, respectively.
SUNSpinWaveTheory.suncouplings
— Methodsuncouplings(J₁₂₃₄::Array{<:Number,4}) -> OperatorSum
Obtain the Couplings for exchange term in Hamiltionian. Like $∑(J_{i1,i2;j3,j4}b†_{i1}b_{i2}b†_{j3}b_{j4})$
SUNSpinWaveTheory.suncouplings
— Methodsuncouplings(Bᵤᵥ::Array{<:Number, 2}) -> MatrixCoupling
Obtain the Couplings for onsite term in Hamiltionian. Like ∑(Bᵢᵤᵥb†ᵢᵤbᵢᵥ )
TightBindingApproximation.commutator
— Methodcommutator(::SUNMagnonic, hilbert::Hilbert{<:Fock{:b}}) -> Diagonal
Get the commutation relation of the Holstein-Primakoff bosons.
SUNSpinWaveTheory.@rotation_gen
— Macro@rotation_gen ::Int -> Vector{Matrix{Complex{Int}}}
Return the generator matrices of SU(n) group.