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

Insert of coupling_rescale in wide_band_bath_discretisation.jl #42

Open
wants to merge 3 commits into
base: main
Choose a base branch
from
Open
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
9 changes: 5 additions & 4 deletions src/diabatic/anderson_holstein.jl
Original file line number Diff line number Diff line change
Expand Up @@ -5,13 +5,14 @@ struct AndersonHolstein{M<:DiabaticModel,B,D,T} <: LargeDiabaticModel
tmp_derivative::Base.RefValue{D}
fermi_level::T
nelectrons::Int
coupling_rescale::Real
end

function AndersonHolstein(model, bath; fermi_level=0.0)
function AndersonHolstein(model, bath; fermi_level=0.0, coupling_rescale=1.0)
tmp_derivative = Ref(NQCModels.zero_derivative(model, zeros(1,1)))
fermi_level = austrip(fermi_level)
nelectrons = count(x -> x <= fermi_level, bath.bathstates)
return AndersonHolstein(model, bath, tmp_derivative, fermi_level, nelectrons)
return AndersonHolstein(model, bath, tmp_derivative, fermi_level, nelectrons, coupling_rescale)
end

NQCModels.nstates(model::AndersonHolstein) = NQCModels.nstates(model.bath) + 1
Expand All @@ -23,7 +24,7 @@ function NQCModels.potential!(model::AndersonHolstein, V::Hermitian, r::Abstract
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)
fillbathcoupling!(V, Vsystem[2,1], model.bath, model.coupling_rescale)
return V
end

Expand All @@ -32,7 +33,7 @@ function NQCModels.derivative!(model::AndersonHolstein, D::AbstractMatrix{<:Herm

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)
fillbathcoupling!(D[I], Dsystem[I][2,1], model.bath, model.coupling_rescale)
end

return D
Expand Down
12 changes: 6 additions & 6 deletions src/diabatic/wide_band_bath_discretisation.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,19 +8,19 @@ function fillbathstates!(out::Hermitian, bath::WideBandBathDiscretisation)
copy!(diagonal, bath.bathstates)
end

function fillbathcoupling!(out::Hermitian, coupling::Real, bath::WideBandBathDiscretisation)
function fillbathcoupling!(out::Hermitian, coupling::Real, bath::WideBandBathDiscretisation, coupling_rescale)
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

do you want to have coupling_rescale in fillbathcoupling! and setcoupling! as an optional keyword? there might be some use cases where these functions might be needed without.

also, we will probably at some point need multiple despatch for these as we will have scenarios where the coupling isn't a global real constant. Therefore, I recommend specifically setting the type to Real already here and in the other similar functions.

first_column = @view out.data[2:end, 1]
setcoupling!(first_column, bath.bathcoupling, coupling)
setcoupling!(first_column, bath.bathcoupling, coupling, coupling_rescale)
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
function setcoupling!(out::AbstractVector, bathcoupling::AbstractVector, coupling::Real, coupling_rescale)
out .= bathcoupling .* coupling .* coupling_rescale
end

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

"""
Expand Down
12 changes: 8 additions & 4 deletions test/wide_band_bath_discretisations.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,10 +8,11 @@ using LinearAlgebra: Hermitian, diagind
bath = NQCModels.TrapezoidalRule(50, -10, 10)
n = NQCModels.nstates(bath)
out = Hermitian(zeros(n+1, n+1))
coupling_rescale = 1.0

NQCModels.DiabaticModels.fillbathstates!(out, bath)
@test all(out[diagind(out)[2:end]] .== bath.bathstates)
NQCModels.DiabaticModels.fillbathcoupling!(out, 1.0, bath)
NQCModels.DiabaticModels.fillbathcoupling!(out, 1.0, bath, coupling_rescale)
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

now you are only testing if the functions run WITH the optional keyword set. you should also have a test where this isn't set.

@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
Expand All @@ -21,27 +22,30 @@ end

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

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

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

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

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

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

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

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