Skip to content

Commit 3b98dfb

Browse files
committed
operator_column for BoseHubbardReal1D2C
1 parent 2f57368 commit 3b98dfb

File tree

3 files changed

+74
-47
lines changed

3 files changed

+74
-47
lines changed

.gitignore

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,3 +1,4 @@
11
/Manifest.toml
22
.DS_Store
33
/test/Manifest.toml
4+
.vscode

src/BoseHubbardReal1D2C.jl

Lines changed: 50 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -48,6 +48,20 @@ function Rimu.starting_address(h::BoseHubbardReal1D2C)
4848
return BoseFS2C(starting_address(h.ha), starting_address(h.hb))
4949
end
5050

51+
struct BoseHubbardReal1D2CColumn{
52+
A,T,O<:BoseHubbardReal1D2C{T}
53+
} <: Rimu.AbstractOperatorColumn{A,T,O}
54+
address::A
55+
hamiltonian::O
56+
diagonal_element::T
57+
end
58+
59+
function Rimu.operator_column(h::BoseHubbardReal1D2C, address::BoseFS2C)
60+
diagonal_element = Rimu.diagonal_element(h, address)
61+
return BoseHubbardReal1D2CColumn(address, h, diagonal_element)
62+
end
63+
64+
5165
Rimu.LOStructure(::Type{<:BoseHubbardReal1D2C{<:Real}}) = IsHermitian()
5266

5367
Base.getproperty(h::BoseHubbardReal1D2C, s::Symbol) = getproperty(h, Val(s))
@@ -56,9 +70,17 @@ Base.getproperty(h::BoseHubbardReal1D2C, ::Val{:hb}) = getfield(h, :hb)
5670
Base.getproperty(h::BoseHubbardReal1D2C{<:Any,<:Any,<:Any,V}, ::Val{:v}) where {V} = V
5771

5872
# number of excitations that can be made
59-
function Rimu.num_offdiagonals(::BoseHubbardReal1D2C, address)
73+
function Rimu.num_offdiagonals(c::BoseHubbardReal1D2CColumn)
74+
address = starting_address(c)
6075
return 2*(num_occupied_modes(address.bsa) + num_occupied_modes(address.bsb))
6176
end
77+
function Rimu.num_offdiagonals(h::BoseHubbardReal1D2C)
78+
address = starting_address(h)
79+
num_particles_a = Rimu.num_particles(address.bsa)
80+
num_particles_b = Rimu.num_particles(address.bsb)
81+
return 2 * (min(num_particles_a, num_modes(h)) +
82+
min(num_particles_b, num_modes(h)))
83+
end
6284

6385
"""
6486
bose_hubbard_2c_interaction(::BoseFS2C)
@@ -79,27 +101,29 @@ end
79101

80102
function Rimu.diagonal_element(ham::BoseHubbardReal1D2C, address::BoseFS2C)
81103
return (
82-
diagonal_element(ham.ha, address.bsa) +
83-
diagonal_element(ham.hb, address.bsb) +
104+
diagonal_element(ham.ha, address.bsa) + # still makes use of the old interface
105+
diagonal_element(ham.hb, address.bsb) + # old interface
84106
ham.v * bose_hubbard_2c_interaction(address)
85107
)
86108
end
87-
88-
function Rimu.get_offdiagonal(ham::BoseHubbardReal1D2C, address, chosen)
89-
nhops = num_offdiagonals(ham,address)
90-
nhops_a = 2 * num_occupied_modes(address.bsa)
91-
if chosen nhops_a
92-
naddress_from_bsa, onproduct = hopnextneighbour(address.bsa, chosen)
93-
elem = - ham.ha.t * onproduct
94-
return BoseFS2C(naddress_from_bsa,address.bsb), elem
95-
else
96-
chosen -= nhops_a
97-
naddress_from_bsb, onproduct = hopnextneighbour(address.bsb, chosen)
98-
elem = -ham.hb.t * onproduct
99-
return BoseFS2C(address.bsa,naddress_from_bsb), elem
100-
end
101-
# return new address and matrix element
102-
end
109+
Rimu.diagonal_element(c::BoseHubbardReal1D2CColumn) = c.diagonal_element
110+
111+
112+
# function Rimu.get_offdiagonal(ham::BoseHubbardReal1D2C, address, chosen)
113+
# nhops = num_offdiagonals(ham,address)
114+
# nhops_a = 2 * num_occupied_modes(address.bsa)
115+
# if chosen ≤ nhops_a
116+
# naddress_from_bsa, onproduct = hopnextneighbour(address.bsa, chosen)
117+
# elem = - ham.ha.t * onproduct
118+
# return BoseFS2C(naddress_from_bsa,address.bsb), elem
119+
# else
120+
# chosen -= nhops_a
121+
# naddress_from_bsb, onproduct = hopnextneighbour(address.bsb, chosen)
122+
# elem = -ham.hb.t * onproduct
123+
# return BoseFS2C(address.bsa,naddress_from_bsb), elem
124+
# end
125+
# # return new address and matrix element
126+
# end
103127

104128
struct OffdiagonalsBoseReal1D2C{
105129
A<:BoseFS2C,T,H<:TwoComponentHamiltonian{T}
@@ -110,9 +134,11 @@ struct OffdiagonalsBoseReal1D2C{
110134
num_hops_a::Int
111135
end
112136

113-
function Rimu.offdiagonals(h::BoseHubbardReal1D2C, a::BoseFS2C)
114-
hops_a = num_offdiagonals(h.ha, a.bsa)
115-
hops_b = num_offdiagonals(h.hb, a.bsb)
137+
function Rimu.offdiagonals(c::BoseHubbardReal1D2CColumn)
138+
h = c.hamiltonian
139+
a = c.address
140+
hops_a = num_offdiagonals(h.ha, a.bsa) # old interface
141+
hops_b = num_offdiagonals(h.hb, a.bsb) # old interface
116142
length = hops_a + hops_b
117143

118144
return OffdiagonalsBoseReal1D2C(h, a, length, hops_a)
@@ -122,10 +148,12 @@ function Base.getindex(s::OffdiagonalsBoseReal1D2C{A,T}, i)::Tuple{A,T} where {A
122148
@boundscheck 1 i s.length || throw(BoundsError(s, i))
123149
if i s.num_hops_a
124150
new_a, matrix_element = get_offdiagonal(s.hamiltonian.ha, s.address.bsa, i)
151+
# old interface
125152
new_add = A(new_a, s.address.bsb)
126153
else
127154
i -= s.num_hops_a
128155
new_b, matrix_element = get_offdiagonal(s.hamiltonian.hb, s.address.bsb, i)
156+
# old interface
129157
new_add = A(s.address.bsa, new_b)
130158
end
131159
return new_add, matrix_element

test/runtests.jl

Lines changed: 23 additions & 25 deletions
Original file line numberDiff line numberDiff line change
@@ -90,12 +90,12 @@ using RimuLegacyHamiltonians: bose_hubbard_2c_interaction
9090

9191
@test hamA == Ĥ2cReal.ha
9292
@test hamB == Ĥ2cReal.hb
93-
@test num_offdiagonals(Ĥ2cReal, aIni2cReal) == 16
94-
@test num_offdiagonals(Ĥ2cReal, aIni2cReal) == num_offdiagonals(Ĥ2cReal.ha, aIni2cReal.bsa) + num_offdiagonals(Ĥ2cReal.hb, aIni2cReal.bsb)
93+
@test num_offdiagonals(operator_column(Ĥ2cReal, aIni2cReal)) == 16
94+
@test 16 == num_offdiagonals(Ĥ2cReal.ha, aIni2cReal.bsa) + num_offdiagonals(Ĥ2cReal.hb, aIni2cReal.bsb)
9595
@test dimension(Ĥ2cReal) == 1225
9696
@test dimension(Float64, Ĥ2cReal) == 1225.0
9797

98-
hp2c = offdiagonals(Ĥ2cReal, aIni2cReal)
98+
hp2c = offdiagonals(operator_column(Ĥ2cReal, aIni2cReal))
9999
@test length(hp2c) == 16
100100
@test hp2c[1][1] == BoseFS2C(BoseFS((0, 2, 1, 1)), BoseFS((1, 1, 1, 1)))
101101
@test hp2c[1][2] -1.4142135623730951
@@ -138,28 +138,26 @@ using RimuLegacyHamiltonians: bose_hubbard_2c_interaction
138138
addr1 = near_uniform(BoseFS2C{1,100,20})
139139
addr2 = near_uniform(BoseFS2C{100,1,20})
140140

141-
for Hamiltonian in (BoseHubbardReal1D2C, BoseHubbardMom1D2C)
142-
@testset "$Hamiltonian" begin
143-
H1 = BoseHubbardReal1D2C(addr1; ta=1.0, tb=2.0, ua=0.5, ub=0.7, v=0.2)
144-
H2 = BoseHubbardReal1D2C(addr2; ta=2.0, tb=1.0, ua=0.7, ub=0.5, v=0.2)
145-
@test starting_address(H1) == addr1
146-
@test LOStructure(H1) == IsHermitian()
147-
148-
hops1 = collect(offdiagonals(H1, addr1))
149-
hops2 = collect(offdiagonals(H2, addr2))
150-
sort!(hops1, by=a -> first(a).bsa)
151-
sort!(hops2, by=a -> first(a).bsb)
152-
153-
addrs1 = first.(hops1)
154-
addrs2 = flip.(first.(hops2))
155-
values1 = last.(hops1)
156-
values2 = last.(hops1)
157-
@test addrs1 == addrs2
158-
@test values1 == values2
159-
160-
@test eval(Meta.parse(repr(H1))) == H1
161-
@test eval(Meta.parse(repr(H2))) == H2
162-
end
141+
@testset "Hamiltonian" begin
142+
H1 = BoseHubbardReal1D2C(addr1; ta=1.0, tb=2.0, ua=0.5, ub=0.7, v=0.2)
143+
H2 = BoseHubbardReal1D2C(addr2; ta=2.0, tb=1.0, ua=0.7, ub=0.5, v=0.2)
144+
@test starting_address(H1) == addr1
145+
@test LOStructure(H1) == IsHermitian()
146+
147+
hops1 = collect(offdiagonals(operator_column(H1, addr1)))
148+
hops2 = collect(offdiagonals(operator_column(H2, addr2)))
149+
sort!(hops1, by=a -> first(a).bsa)
150+
sort!(hops2, by=a -> first(a).bsb)
151+
152+
addrs1 = first.(hops1)
153+
addrs2 = flip.(first.(hops2))
154+
values1 = last.(hops1)
155+
values2 = last.(hops1)
156+
@test addrs1 == addrs2
157+
@test values1 == values2
158+
159+
@test eval(Meta.parse(repr(H1))) == H1
160+
@test eval(Meta.parse(repr(H2))) == H2
163161
end
164162
end
165163
@testset "1D Bosons (2-component)" begin

0 commit comments

Comments
 (0)