CPIMC
Base.:+ — MethodBase.:+(x::UpdateCounter, y::UpdateCounter)Addition of counters for example to add counters from different Markov-chains.
Base.getindex — MethodBase.getindex(kinks::Kinks{T}, τ::ImgTime) where {T}Allow the kinks[τ]-syntax to retrieve a Kink at a specific ImgTime.
Base.haskey — MethodBase.haskey(ck::Kinks, key::ImgTime)Check if a Kinks-object contains a kink at a specific time.
CPIMC.Woffdiag_element — MethodWoffdiag_element(::Model, ::Ensemble, ::Orbital, ::Orbital, ::Orbital, ::Orbital)Return the offdiagonal many body matrix element of the interaction operator
\[\frac{1}{4} ( w_{ijkl} - w_{ijlk} ) (\pm \langle\{\tilde{n}\}|\{n\}^{ij}_{kl}\rangle)\]
for an excitation given by creating orbitals i,j and annihilating orbitals k, l as given by the Slater-Condon rules.
CPIMC.add_kinks! — Methodadd_kinks!(ck::Kinks, s)
Add the kinks in s to ck.
CPIMC.add_kinks — Methodadd_kinks(ck::Kinks, s)Return a copy of the kinks in ck with the kinks in s added.
CPIMC.add_orbs! — Methodadd_orbs!(occ, s)Add the orbitals in s to the occupation occ.
CPIMC.add_orbs — Methodadd_orbs(occ, orbs)Return a copy of the occupation occ with the orbitals in orbs added.
CPIMC.adjacent_kinks — Methodadjacent_kinks(itr::Any, ::Any)Return a tuple of the kink in itr that is closest right of the time in the second argument. This expects that methods for the functions next() and prev() are defined for the given argument types. Used for getting a tuple of the neighbouring kinks from some imaginary time.
CPIMC.adjacent_kinks_affecting_orbs — Methodadjacent_kinks_affecting_orbs(::Any, ::Any, ::Any)Return a tuple of the closest kink to the right and the closest kink to the left of some imaginary time that affect one of the orbitals in some collection. If there are no such kinks, return (nothing,nothing). The first argument is expected to be an iterable which contains pairs as elements, with the first element of each pair containing an imaginary time and the second element of each pair containing a kink. The second argument is expected to be a collection of orbitals, that are to be considered regarding affection with one of the kinks in the first argument. The third argument is expected to be an imaginary time from which the neighbouring kinks affected by any orbitals from the second argument are to be determined.
CPIMC.annihilators — Methodannihilators(x::T4)return a set of the two annihilators which are affected by a T4 kink
CPIMC.apply_step! — Methodapply_step!(c::Configuration, step::Step)
apply_step!(c::Configuration, steps)Apply the changes given by the second argument to a Configuration.
CPIMC.apply_step — Methodapply_step(c::Configuration, step::Step)
apply_step(c::Configuration, steps)Return the result of applying the changes given by the second argument to a Configuration.
CPIMC.basis — Methodbasis(c::Configuration{T})Return type parameter T of the configuration.
CPIMC.creators — Methodcreators(::T4)return a set of the two creators which are affected by a T4 kink
CPIMC.drop_kinks! — Methoddrop_kinks!(ck, s)
Remove the kinks in s from ck.
CPIMC.drop_kinks — Methoddrop_kinks(ck, s)Return a copy of the kinks in ck without the kinks in s.
CPIMC.drop_orbs! — Methoddrop_orbs!(occ, s)Remove the orbitals in s from the occupation occ.
CPIMC.drop_orbs — Methoddrop_orbs(occ, orbs)Return a copy of the occupation occ without the orbitals in orbs.
CPIMC.excitations — Methodexcitations(ck::Kinks)Return a list of the excitations of a Kinks-object, used to allow the use of dictionary syntax.
CPIMC.excite! — Methodexcite!(::Set{T}, i::T, j::T, k::T, l::T) where {T}
excite!(o::Set{T}, κ::T4{T})Apply an excitation in-place to a set of basis states that is given by creating orbitals i, j and annihilating the orbitals k, l.
CPIMC.excite — Methodexcite(o::Set{T}, κ::Pair{ImgTime,T4{T}})Apply a T4-kink to a set of basis states for a pair of a time and a T4-kink. This is useful for iterating over a Kinks-object when calculating the occupations at a specific time.
CPIMC.excite — Methodexcite(::Set{T}, ::T4{T})Apply a kink to a set of basis states, i.e. return a set of basis states where the states specified by the creators of the kink are added and the states specified by the annihilators of the kink are dropped from the given set of basis states. This has the physical meaning of a two-particle scattering event where two (quasi-)particles in a many-body state change the single-particle states they occupy. In occupation number representation this reads
\[a^{\dagger}_i a^{\dagger}_j a_k a_l |\{n\}\rangle\]
for creator orbitals i and j and annihilator orbitals k and l. This function assumes fermionic particle statistics (Pauli principle),
\[a^{\dagger}_i a^{\dagger}_i = a_i a_i = 0\]
i.e. the target (creator) states must not be occupied and the initial (annihilator) states must be occupied in the given set of states.
CPIMC.in_open_interval — Methodin_open_interval(τ, τ_first, τ_last)Check if τ lies in the open ImgTime-interval (τ_first,τlast), assuming thatτfirstis the left border andτ_last` the right one.
CPIMC.is_type_1 — Methodis_type_1(left_kink::T4, right_kink::T4)Returns True if leftkink and rightkink are entangled in a Type-1 way. This does not check wether the two kinks are neighbouring.
CPIMC.isunaffected — Methodisunaffected(Kinks, orb)Return if an orbital is not affected by any kink.
CPIMC.isunaffected_in_interval — Methodisunaffected_in_interval(kinks, orb, τ_first::ImgTime, τ_last::ImgTime)Return if an orbital is not affected by any of the kinks from ck in the open interval (τ_first,τ_last).
CPIMC.kinks_affecting_orbs — Method kinks_affecting_orbs(itr::Any, itr::Any)Return all kinks in the first argument that affect one ore more of the orbitals the second argument. The first itr is expected to contain tuples or at least types on which last() is defined and last(i) is intended to contain a ::T4 for all elements i ∈ itr or at least a type which has fields with the names i, j, k, l.
CPIMC.kinks_from_periodic_interval — Methodkinks_from_periodic_interval(::Kinks, τ1, τ2)return kinks with τ ∈ (τ1,τ2) if τ1 < τ2 and τ ∈ (τ2,1) ∪ (0,τ1) if τ1 > τ2
CPIMC.ladder_operator_order_factor — Methodladder_operator_order_factor(sortedlist::Array{<:Orbital,1})Returns $1$ or $-1$ depending on the order of all ladder operators as given by a list of orbitals.
CPIMC.ladder_operator_order_factor — Methodladder_operator_order_factor(ck::Kinks)Returns $1$ or $-1$ depending on the order of all ladder operators as given by the orbitals that affect each kink in the time-ordering of the kinks and in the conventional ordering i, j, k, l, used in the sign estimator.
CPIMC.left_type_1_chain_length — Functionleft_type_1_chain_length(ck::Kinks, τ, counted_τs = [])Returns the length of the chain of type-1-entaglements starting with the Kink at τ counting to the left.
CPIMC.longest_type_1_chain_length — Methodlongest_type_1_chain_length(ck::Kinks)Returns the longest chain of type-1-entaglements in ck.
CPIMC.measure! — Methodmeasure!(m::Model, e::Ensemble, c::Configuration, estimators)Perform measurements on a Configuration. estimators needs to be a dictionary that contains tuples (::OnlineStat, ::Function).
CPIMC.next — Methodnext(a, τ)Return index of the first kink after to the given τ. If there is no kink with ::ImgTime larger than τ, return index of first kink. Throw a DomainError if an empty list is passed as first argument.
This function assumes that the kinks in a are ordered with respect to their imaginary time!
CPIMC.next_affecting — Method next_affecting(a::Any, orbs::Any, τ::Any)Return the index of the closest kink to the right of τ that affects one of the orbitals in orbs. If no orbital in orbs is affected by any kink return 0.
CPIMC.occupations_at — Methodoccupations_at(c::Configuration, τ::ImgTime)Return the occupied orbitals to the right of τ.
CPIMC.occupations_at — Methodoccupations_at(o::Set{T}, kinks::Kinks{T})Return the occupied orbitals after applying all kinks to initial occupation.
CPIMC.orbs — Methodorbs(::T2)Return a Tuple of all orbitals that are affected by a T2-kink in the conventional ordering i, j.
CPIMC.orbs — Methodorbs(::T4)Return a Tuple of all orbitals that are affected by a T4-kink in the conventional ordering i, j, k, l. This corresponds to the matrix element $w_{ijkl}$ in the same order.
CPIMC.prev — Methodprev(a, τ)Return the index of the first kink earlier than τ::ImgTime. If there is no such kink the last kink is returned. Throw a DomainError if an empty list is passed as first argument.
This function assumes that the kinks in x are ordered with respect to their imaginary time!
CPIMC.prev_affecting — Method prev_affecting(a::Any, orbs::Any, τ::Any)Return the index of the closest kink to the left of τ that affects one of the orbitals in orbs. If no orbital in orbs is affected by any kink return 0.
CPIMC.print_rates — Methodprint_rates(dict::Dict{Any,UpdateCounter})Calculate and print acceptance statistics.
CPIMC.print_results — Methodprint_results(measurements, e::Ensemble)Print all measured observables, potentially taking the sign into account.
CPIMC.random_shuffle — Methodrandom_shuffle(::Any, ::Any)shuffle or not shuffle the given tuple with equal propability: Return either the same tuple, or the tuple in reverse order, with equal probability 0.5.
CPIMC.right_type_1_chain_length — Functionright_type_1_chain_length(ck, τ, count = 0)Returns the length of the chain of type-1-entaglements starting with the Kink at τ counting to the right.
CPIMC.right_type_1_count — Methodright_type_1_count(ck::Kinks)Returns the number of kinks that are type 1 entagled with their right neighbour.
CPIMC.shuffle_annihilators — Methodshuffle_annihilators(::Any, ::Any, ::Any, ::Any)Return a tuple where the last two arguments are randomly shuffled with equal probability 0.5 in comparison to the ordering of the input.
CPIMC.shuffle_creators — Methodshuffle_creators(::Any, ::Any, ::Any, ::Any)Return a tuple where the first two arguments are randomly shuffled with equal probability 0.5 in comparison to the ordering of the input.
CPIMC.shuffle_indices — Methodshuffle_indices(::Any, ::Any, ::Any, ::Any)Return a tuple where both the first two arguments and the last two arguments are randomly shuffled with equal probability 0.5 respectively in comparison to the ordering of the input.
CPIMC.sign_offdiagonal_product — Methodsign_offdiagonal_product(::Model, ::Configuration)Return the sign of the product of two-particle terms in the offdiagonal many-body matrix elements. These are given by function wminus, i.e.
\[\langle \{\tilde{n}\} | a^{\dagger}_i a^{\dagger}_j a_k a_l | \{n\} \rangle = \pm ( w_{ijkl} - w_{ijlk} ) \text{ for } \{\tilde{n}\} = \{n\}_{kl}^{ij}\]
the term in the braces may be negative and this function returns the product of the sign of all these contributions $w_{ijkl} - w_{ijlk}$ from all kinks. The sign $\pm$ is determined from the permutation factor of the orbitals i,j,k,l and is not calculated here.
used for the calculation the sign of the weight function
CPIMC.signum — Methodsignum(m::Model, c::Configuration)Calculate the sign of the Configuration's weight.
CPIMC.sweep! — Methodsweep!(m::Model, e::Ensemble, c::Configuration, updates, estimators, steps::Int, sampleEvery::Int, throwAway::Int)Generate a markov chain of length steps using the Metropolis-Hastings algorithm with the updates given in updates. After throwAway steps have been performed, the observables given in estimators are calculated every sampleEvery steps.
Returns a dictionary containing an UpdateCounter for every function given in updates.
CPIMC.time_ordered_orbs — Methodtime_ordered_orbs(ck::Kinks{T}) where {T <: Orbital}Get a list of orbitals that affect each kink in the time-ordering of the kinks and in the conventional ordering i, j, k, l.
CPIMC.times — Methodtimes(ck::Kinks)Return a list of the imaginary-times of a Kinks-object, used to allow the use of dictionary syntax.
CPIMC.times_from_periodic_interval — Methodtimes_from_periodic_interval(::Kinks, ::ImgTime, ::ImgTime)return a list of all times of kinks with τ ∈ (τ1,τ2) if τ1 < τ2 or τ ∈ (τ2,1) ∪ (0,τ1) if τ1 > τ2 in the periodic ordering suggested by the relation of the first time-argument τ1 to the second time-argument τ2
CPIMC.update! — Methodupdate!(m::Model, e::Ensemble, c::Configuration, updates)Obtain the next state of the Markov chain by randomly selecting an element of updates and applying the proposed changes to c with the calculated probability. Returns a tuple containing the index of the chosen function and :accept, :reject or :trivial if the function yielded an empty Step.
CPIMC.wminus — Methodwminus(::Model, i,j,k,l)return the difference of the two-particle matrix element with the same but the last two indices transposed antisymmetric difference of the two-particle matrix elements:
\[w^-_{ijkl} = w_{ijkl} - w_{ijlk}\]
This is called the antisymmetrized two-particle matrix element. This is an abbreviation for these terms arising in the Slater-Condon rules for the calculation of the many-body matrix elements via the one- and two-particle matrix elements of the underlying single-particle basis.
CPIMC.Δ — MethodΔ(::ImgTime,::ImgTime)return the length of the periodic interval between the two ::ImgTimes
the periodic distance of two imaginary times:
`Δ(τ1,τ2) = τ2 - τ1` if `τ2 >= τ1` and
`Δ(τ1,τ2) = 1 - (τ1 - τ2) = 1 + τ2 - τ1` elseCPIMC.ΔT_element — MethodΔT_element(::Model, i,j,k,l)return the change in the kinetic many body matrix element due to creating orbitals i, j and annihilating orbitals k, l
CPIMC.ΔW_diag — MethodΔW_diag(m::Model, i, j, k, l, occ)change in the diagonal interaction matrix element due to a change in the occupation occ in a periodic interval (τ1,τ2) where no kinks occur the change in the occupation is assumed to consist in a creation of two orbitals i, j and in the annihlation of two orbitals k, l
CPIMC.ΔW_diag — MethodΔW_diag(m::Model, i, j, occ)change in the diagonal interaction matrix element due to a change in the occupation occ in a periodic interval (τ1,τ2) where no kinks occur the change in the occupation is assumed to consist in a creation of one orbitals i and in the annihlation of one orbitals j
CPIMC.ΔWdiag_element — MethodΔWdiag_element(::Model, ::Ensemble, ::Configuration, i, j, k, l, τ1, τ2)Calculate the change in the diagonal interaction many-body matrix element due to a change in the occupations given by creating two orbitals i and j and annihilating two orbitals k, l in the interval (τ1, τ2). This interval may be periodically extended over the bounds (0,1) if τ1 > τ2, i.e. the change in the occupation is considered for (τ1,1] ∪ [0,τ2) in that case. We do not need to evaluate the diagonal interaction between all orbitals in all time-intervalls, but it is sufficient to evaluate the full diagonal interaction with the occupations at the start of the intervall and then consider only contributions of orbitals that are changed by kinks in the intervall.
CPIMC.ΔWdiag_element — MethodΔWdiag_element(::Model, ::Ensemble, ::Configuration, i, j, τ1, τ2)Calculate the change in the diagonal interaction many-body matrix element due to a change in the occupations given by creating an orbital i and annihilating an orbital j in the interval (τ1, τ2). This interval may be periodically extended over the bounds (0,1) if τ1 > τ2, i.e. the change in the occupation is considered for (τ1,1] ∪ [0,τ2) in that case. We do not need to evaluate the diagonal interaction between all orbitals in all time-intervalls, but it is sufficient to evaluate the full diagonal interaction with the occupations at the start of the intervall and then consider only contributions of orbitals that are changed by kinks in the intervall.
CPIMC.ΔWoffdiag_element — MethodΔWoffdiag_element(::Model, e::Ensemble, itr, itr)return the change in the offdiagonal many body matrix element of the interaction operator given by adding the kinks in the first iterable and removing the kinks in the second iterable
both arguments itr are required to be iterables containing kinks
CPIMC.τ_borders — Method τ_borders(::Kinks{T}, orbs, ::ImgTime) where {T <: Orbital}Return a tuple of the ImgTime of the closest kink to the right and the ImgTime of the closest kink to the left of τ that affect one of the orbitals in orbs. If no orbital in os is affected by and kink from the collection in the first argument, return a tuple of the interval bounds (ImgTime(0), ImgTime(1)).
CPIMC.τ_next_affecting — Method τ_next_affecting(::Any, ::Any, ::Any)Return the ImgTime of the closest kink to the right of τ that affects one of the orbitals in os.
CPIMC.τ_prev_affecting — Method τ_prev_affecting(::Any, ::Any, ::Any)Return the ImgTime of the closest kink to the left of τ that affects one of the orbitals in os.
Estimators
CPIMC.Estimators.E — MethodE(m::Model, e::Ensemble, c::Configuration)estimator for the interaction energy This estimator is redundant because we can calculate the full energy as the sum of the kinetic and the interaction part.
CPIMC.Estimators.Ekin — MethodEkin(m::Model, e::Ensemble, c::Configuration)estimator for the kinetic energy
CPIMC.Estimators.K — MethodK(m::Model, e::Ensemble, c::Configuration)estimator for the number of kinks
CPIMC.Estimators.W — MethodW(m::Model, e::Ensemble, c::Configuration)estimator for the interaction energy This estimator is redundant because we can calculate the interaction energy from Kfermion and Wdiag.
CPIMC.Estimators.W_diag — MethodW_diag(m::Model, e::Ensemble, c::Configuration)estimator for the diagonal contribution to the interaction energy
CPIMC.Estimators.W_off_diag — MethodW_off_diag(m::Model, e::Ensemble, c::Configuration)estimator for the offdiagonal contribution to the interaction energy This estimator is redundant because we can calculate the offdiagonal interaction energy from K_fermion.
CPIMC.Estimators.occupations — Functionoccupations(m::Model, e::Ensemble, c::Configuration, emax::Int=100) :: Array{Float64,1}estimator for the occupations of the emax:Int lowest single particle energy eigenvalues
CPIMC.signum — Methodsignum(m::Model, e::Ensemble, c::Configuration)estimator for the sign of the weight function
PlaneWaves
CPIMC.PlaneWaves.dimension — Methoddimension(o::PlaneWave{D}) where {D}return the dimension of the momentum vector of an Orbital
CPIMC.PlaneWaves.dimension — Methoddimension(os::Set{<:PlaneWave{D}}) where {D}return the dimension of the momentum vectors of a Set{<:Orbital{D}}
CPIMC.PlaneWaves.find_fourth_orb_for_kink — Methodfind_fourth_orb_for_kink(same_kind_ladder_operator, other_kind_ladder_operator1, other_kind_ladder_operator2)Find the fourth orb to build a Kink accoring to spin/momentum conservation, where the Operator thats acts on the fourth orb is of the same kind (creator/ annihilator) as the one that acts on samekindladderoperator, While on otherkindladderoperator1 and otherkindladder_operator2 the other kind of operators will act.
CPIMC.PlaneWaves.flip — Methodflip(o::PlaneWave)return an PlaneWave with the same momentum vector but opposite spin projection as the given o::PlaneWave
CPIMC.PlaneWaves.flip — Methodflip(s::Spin)return opposite spin projection
CPIMC.PlaneWaves.fractional_spin_polarization — Methodfractional_spin_polarization(c::Configuration{<:PlaneWave})calculate the fractional spin polarization from a set of Orbitals which each must have a field spin of type @enum Spin Down Up, the fractional spin polarization is defined as the ratio of the absolute difference of the number of particles which Spin Down and Spin Up and the particle number, thus here no convention is made as to which spin component should be occupied more often, as opposed to some literature where N↑ > N↓ is used
CPIMC.PlaneWaves.fractional_spin_polarization — Methodfractional_spin_polarization(occ::Set{PlaneWave)calculate the fractional spin polarization from a set of Orbitals which each must have a field spin of type @enum Spin Down Up, the fractional spin polarization is defined as the ratio of the absolute difference of the number of particles which Spin Down and Spin Up and the particle number, thus here no convention is made as to which spin component should be occupied more often, as opposed to some literature where N↑ > N↓ is used
CPIMC.PlaneWaves.sphere — Methodsphere(o::PlaneWave{1}; dk=2)
sphere(o::PlaneWave{2}; dk=2)
sphere(o::PlaneWave{3}; dk=2)Return the set of all PlaneWaves in a sphere of radius dk around o in momentum space.
CPIMC.PlaneWaves.sphere_with_same_spin — Methodsphere_with_same_spin(o::PlaneWave{1}; dk=2)
sphere_with_same_spin(o::PlaneWave{2}; dk=2)
sphere_with_same_spin(o::PlaneWave{3}; dk=2)Return the set of all PlaneWaves in a sphere of radius dk around o in momentum space that also have the same spin as o.
CPIMC.PlaneWaves.ξ — Methodξ(c::Configuration{<:PlaneWave})calculate the fractional spin polarization from a set of Orbitals which each must have a field spin of type @enum Spin Down Up, the fractional spin polarization is defined as the ratio of the absolute difference of the number of particles which Spin Down and Spin Up and the particle number, thus here no convention is made as to which spin component should be occupied more often, as opposed to some literature where N↑ > N↓ is used
CPIMC.PlaneWaves.ξ — Methodξ(occ::Set{PlaneWave})calculate the fractional spin polarization from a set of Orbitals which each must have a field spin of type @enum Spin Down Up, the fractional spin polarization is defined as the ratio of the absolute difference of the number of particles which Spin Down and Spin Up and the particle number, thus here no convention is made as to which spin component should be occupied more often, as opposed to some literature where N↑ > N↓ is used
UniformElectronGas
CPIMC.UniformElectronGas.EF — MethodEF(rs, ξ, d::Int)calculate the Fermi energy in Hartree a.u. from Brueckner parameter rs for the density, the fractional spin polarization ξ and spatial dimension d
CPIMC.UniformElectronGas.internal_energy_factor — Methodinternal_energy_factor(N, rs, d)This factor is used to convert an energy quantity from Hartree a.u. to internal units. A value of the quantity in Hartree a.u. has to be multiplied by this factor to yield the corresponding value of the quantity in internal units. A value of the quantity in internal units has to be divided by this factor to yield the corresponding value of the quantity in Hartree a.u.
CPIMC.UniformElectronGas.kF — MethodkF(rs, ξ, d::Int)calculate the Fermi wavenumber in Hartree a.u. from Brueckner parameter rs for the density, the fractional spin polarization ξ and spatial dimension d
CPIMC.UniformElectronGas.rs — Methodrs(N, λ, d::Int)calculate the Brueckner parameter rs for the density from the particle number N and λ, which is the coupling constant of the Hamiltonian in internal units
CPIMC.UniformElectronGas.β — Methodβ(Θ, rs, N, ξ, d::Int)calculate β in internal units
CPIMC.UniformElectronGas.β_Ha — Methodβ_Ha(Θ, EF)calculate the inverse temperature in Hartree a.u. from the reduced temperature Θ and the Fermi energy EF in Hartree a.u. the inverse temperature is defined by β = 1 / (kB*T) where kB is the Boltzmann constant and T is the temperature
CPIMC.UniformElectronGas.λ — Methodλ(N, rs, d::Int)calculate λ, which is the coupling constant of the Hamiltonian in internal units from the particle number N, Bruecker parameter rs for the density and spatial dimension d the differences for d ∈ {2,3} reflect the different expressions for the Coulomb interaction, whichs Fourier component is vq = rac{4π}{q^2} in 3D and vq = rac{2π}{q} in 2D the one-dimensional case involves a short-distance cutoff to avoid the non-integrable divergence of the interaction at zero distance and is not implemented
CPIMC.energy — Methodenergy(m::UEG, o::PlaneWave)The single-particle energy of a plane wave with momentum $\vec{k}$ is given by (Hartree atomic units):
$ \epsilon_k = \vec{k}^2 / 2. $
In the internal unit system we use integer $k$-values, which results in:
$ \epsilon_k = \vec{k}^2. $
CPIMC.w — Methodw(i::Orbital, j::Orbital, k::Orbital, l::Orbital)Return the two-particle matrix element of the Coulomb interaction. Zero is returned if momentum or spin is not conserved, in consequence of the Bloch-theorem and the spin kronecker-delta in the plane-spin-wave basis.
The singular contribution (i.vec == k.vec) is also removed, i.e. set to zero. This property is not a result of the Bloch theorem for the two-particle Coulomb matrix element in the plane wave basis but is added here to allow for straightforward summation over all combination of orbitals without explictly excluding this component.
It is not part of the UEG many-body Hamiltonian due to the physical assumption that contributions from the uniform background cancel with this diverging term in the thermodynamic limit.