Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add AndersonHolstein model #19

Merged
merged 3 commits into from
Jun 29, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 2 additions & 0 deletions Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,7 @@ authors = ["James Gardner <james.gardner1421@gmail.com>"]
version = "0.8.10"

[deps]
FastGaussQuadrature = "442a2c76-b920-505d-bb47-c5924d526838"
ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210"
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
NQCBase = "78c76ebc-5665-4934-b512-82d81b5cbfb7"
Expand All @@ -18,6 +19,7 @@ UnitfulAtomic = "a7773ee8-282e-5fa2-be4e-bd808c38a91a"
UnitfulRecipes = "42071c24-d89e-48dd-8a24-8a12d9b8861f"

[compat]
FastGaussQuadrature = "0.4"
ForwardDiff = "0.10"
NQCBase = "0.1, 0.2"
Parameters = "0.12"
Expand Down
9 changes: 9 additions & 0 deletions src/diabatic/DiabaticModels.jl
Original file line number Diff line number Diff line change
Expand Up @@ -161,4 +161,13 @@ export ErpenbeckThoss
include("widebandbath.jl")
export WideBandBath

include("wide_band_bath_discretisation.jl")
export TrapezoidalRule
export ShenviGaussLegendre
export ReferenceGaussLegendre
export FullGaussLegendre

include("anderson_holstein.jl")
export AndersonHolstein

end # module
43 changes: 43 additions & 0 deletions src/diabatic/anderson_holstein.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,43 @@

struct AndersonHolstein{M<:DiabaticModel,B} <: LargeDiabaticModel
model::M
bath::B
end

NQCModels.nstates(model::AndersonHolstein) = NQCModels.nstates(model.bath) + 1
NQCModels.ndofs(model::AndersonHolstein) = NQCModels.ndofs(model.model)
NQCModels.nelectrons(model::AndersonHolstein) = fld(NQCModels.nstates(model.bath), 2)
NQCModels.fermilevel(::AndersonHolstein) = 0

function NQCModels.potential!(model::AndersonHolstein, V::Hermitian, r::AbstractMatrix)
Vsystem = NQCModels.potential(model.model, r)
V[1,1] = Vsystem[2,2] - Vsystem[1,1]
fillbathstates!(V, model.bath)
fillbathcoupling!(V, Vsystem[2,1], model.bath)
return V
end

function NQCModels.derivative!(model::AndersonHolstein, D::AbstractMatrix{<:Hermitian}, r::AbstractMatrix)

Dsystem = NQCModels.derivative(model.model, r)

for I in eachindex(Dsystem, D)
D[I][1,1] = Dsystem[I][2,2] - Dsystem[I][1,1]
fillbathcoupling!(D[I], Dsystem[I][2,1], model.bath)
end

return D
end

function NQCModels.state_independent_potential(model::AndersonHolstein, r::AbstractMatrix)
Vsystem = NQCModels.potential(model.model, r)
return Vsystem[1,1]
end

function NQCModels.state_independent_derivative!(model::AndersonHolstein, ∂V::AbstractMatrix, r::AbstractMatrix)
Dsystem = NQCModels.derivative(model.model, r)
for I in eachindex(∂V, Dsystem)
∂V[I] = Dsystem[I][1,1]
end
return ∂V
end
138 changes: 138 additions & 0 deletions src/diabatic/wide_band_bath_discretisation.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,138 @@
using FastGaussQuadrature: gausslegendre

abstract type WideBandBathDiscretisation end
NQCModels.nstates(bath::WideBandBathDiscretisation) = length(bath.bathstates)

function fillbathstates!(out::Hermitian, bath::WideBandBathDiscretisation)
diagonal = view(out, diagind(out)[2:end])
copy!(diagonal, bath.bathstates)
end

function fillbathcoupling!(out::Hermitian, coupling::Real, bath::WideBandBathDiscretisation)
first_column = @view out.data[2:end, 1]
setcoupling!(first_column, bath.bathcoupling, coupling)
first_row = @view out.data[1, 2:end]
copy!(first_row, first_column)
end

function setcoupling!(out::AbstractVector, bathcoupling::AbstractVector, coupling::Real)
out .= bathcoupling .* coupling
end

function setcoupling!(out::AbstractVector, bathcoupling::Real, coupling::Real)
fill!(out, bathcoupling * coupling)
end

"""
TrapezoidalRule{B,T} <: WideBandBathDiscretisation

Discretise wide band continuum using trapezoidal rule.
Leads to evenly spaced states and constant coupling.
"""
struct TrapezoidalRule{B,T} <: WideBandBathDiscretisation
bathstates::B
bathcoupling::T
end

function TrapezoidalRule(M, bandmin, bandmax)
ΔE = bandmax - bandmin
bathstates = range(bandmin, bandmax, length=M)
coupling = sqrt(ΔE / M)
return TrapezoidalRule(bathstates, coupling)
end

"""
ShenviGaussLegendre{T}

Defined as described by Shenvi et al. in J. Chem. Phys. 130, 174107 (2009).
The position of the negative sign for the state energy level has been moved to ensure the states are sorted from lowest to highest.
"""
struct ShenviGaussLegendre{T} <: WideBandBathDiscretisation
bathstates::Vector{T}
bathcoupling::Vector{T}
end

function ShenviGaussLegendre(M, bandmin, bandmax)
M % 2 == 0 || throw(error("The number of states `M` must be even."))
knots, weights = gausslegendre(div(M, 2))
ΔE = bandmax - bandmin

bathstates = zeros(M)
for i in eachindex(knots)
bathstates[i] = ΔE/2 * (-1/2 + knots[i]/2)
end
for i in eachindex(knots)
bathstates[i+length(knots)] = ΔE/2 * (1/2 + knots[i]/2)
end

bathcoupling = zeros(M)
for (i, w) in zip(eachindex(bathcoupling), Iterators.cycle(weights))
bathcoupling[i] = sqrt(ΔE * w) / 2
end

return ShenviGaussLegendre(bathstates, bathcoupling)
end

"""
ReferenceGaussLegendre{T}

Implementation translated from Fortran code used for simulations of Shenvi et al. in J. Chem. Phys. 130, 174107 (2009).
Two differences from ShenviGaussLegendre:
- Position of minus sign in energy levels has been corrected.
- Division by sqrt(ΔE) in the coupling.
"""
struct ReferenceGaussLegendre{T} <: WideBandBathDiscretisation
bathstates::Vector{T}
bathcoupling::Vector{T}
end

function ReferenceGaussLegendre(M, bandmin, bandmax)
M % 2 == 0 || throw(error("The number of states `M` must be even."))
@warn "(ReferenceGaussLegendre) The division by sqrt(ΔE) in the coupling is incorrect. Use ShenviGaussLegendre instead."
knots, weights = gausslegendre(div(M, 2))
ΔE = bandmax - bandmin

bathstates = zeros(M)
for i in eachindex(knots)
bathstates[i] = ΔE/2 * (-1/2 + knots[i]/2)
end
for i in eachindex(knots)
bathstates[i+length(knots)] = ΔE/2 * (1/2 + knots[i]/2)
end

bathcoupling = zeros(M)
for (i, w) in zip(eachindex(bathcoupling), Iterators.cycle(weights))
bathcoupling[i] = sqrt(w) / 2
end

return ReferenceGaussLegendre(bathstates, bathcoupling)
end

"""
FullGaussLegendre{T} <: WideBandBathDiscretisation

Use Gauss-Legendre quadrature to discretise the continuum across the entire band width.
This is similar to the ShenviGaussLegendre except that splits the continuum at the Fermi level into two halves.
"""
struct FullGaussLegendre{T} <: WideBandBathDiscretisation
bathstates::Vector{T}
bathcoupling::Vector{T}
end

function FullGaussLegendre(M, bandmin, bandmax)

knots, weights = gausslegendre(M)
ΔE = bandmax - bandmin

bathstates = zeros(M)
for i in eachindex(knots, bathstates)
bathstates[i] = ΔE/2 * knots[i]
end

bathcoupling = zeros(M)
for i in eachindex(bathcoupling, weights)
bathcoupling[i] = sqrt(weights[i] * ΔE/2)
end

return FullGaussLegendre(bathstates, bathcoupling)
end
6 changes: 6 additions & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,9 @@ using NQCBase
using NQCModels
using LinearAlgebra
using FiniteDiff
using SafeTestsets

@time @safetestset "Wide band bath discretisations" begin include("wide_band_bath_discretisations.jl") end

function finite_difference_gradient(model::NQCModels.AdiabaticModels.AdiabaticModel, R)
f(x) = potential(model, x)
Expand Down Expand Up @@ -91,6 +94,9 @@ end
@test test_model(ErpenbeckThoss(Γ=2.0), 1)
@test test_model(WideBandBath(ErpenbeckThoss(Γ=2.0); step=0.1, bandmin=-1.0, bandmax=1.0), 1)
@test test_model(WideBandBath(GatesHollowayElbow(); step=0.1, bandmin=-1.0, bandmax=1.0), 2)
@test test_model(AndersonHolstein(ErpenbeckThoss(Γ=2.0), TrapezoidalRule(10, -1, 1)), 1)
@test test_model(AndersonHolstein(ErpenbeckThoss(Γ=2.0), ShenviGaussLegendre(10, -1, 1)), 1)
@test test_model(AndersonHolstein(GatesHollowayElbow(), ShenviGaussLegendre(10, -1, 1)), 2)
end

@testset "FrictionModels" begin
Expand Down
47 changes: 47 additions & 0 deletions test/wide_band_bath_discretisations.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,47 @@
using Test
using NQCModels

using Unitful
using LinearAlgebra: Hermitian, diagind

@testset "TrapezoidalRule" begin
bath = NQCModels.TrapezoidalRule(50, -10, 10)
n = NQCModels.nstates(bath)
out = Hermitian(zeros(n+1, n+1))

NQCModels.DiabaticModels.fillbathstates!(out, bath)
@test all(out[diagind(out)[2:end]] .== bath.bathstates)
NQCModels.DiabaticModels.fillbathcoupling!(out, 1.0, bath)
@test all(isapprox.(out[1,2:end], 1.0 / sqrt(50/20)))
@test all(isapprox.(out[2:end,1], 1.0 / sqrt(50/20)))
end

@testset "ShenviGaussLegendre" begin
bath = NQCModels.ShenviGaussLegendre(20, -10, 10.0)

n = NQCModels.nstates(bath)
out = Hermitian(zeros(n+1, n+1))

NQCModels.DiabaticModels.fillbathstates!(out, bath)
NQCModels.DiabaticModels.fillbathcoupling!(out, 1.0, bath)
end

@testset "ReferenceGaussLegendre" begin
bath = NQCModels.ReferenceGaussLegendre(20, -10, 10.0)

n = NQCModels.nstates(bath)
out = Hermitian(zeros(n+1, n+1))

NQCModels.DiabaticModels.fillbathstates!(out, bath)
NQCModels.DiabaticModels.fillbathcoupling!(out, 1.0, bath)
end

@testset "FullGaussLegendre" begin
bath = NQCModels.FullGaussLegendre(20, -10, 10.0)

n = NQCModels.nstates(bath)
out = Hermitian(zeros(n+1, n+1))

NQCModels.DiabaticModels.fillbathstates!(out, bath)
NQCModels.DiabaticModels.fillbathcoupling!(out, 1.0, bath)
end