We compute the polarizability of a Helium atom. The polarizability is defined as the change in dipole moment $$ μ = ∫ r ρ(r) dr $$ with respect to a small uniform electric field $E = -x$.
We compute this in two ways: first by finite differences (applying a finite electric field), then by linear response. Note that DFTK is not really adapted to isolated atoms because it uses periodic boundary conditions. Nevertheless we can simply embed the Helium atom in a large enough box (although this is computationally wasteful).
As in other tests, this is not fully converged, convergence parameters were simply selected for fast execution on CI,
using DFTK
using LinearAlgebra
a = 10.
lattice = a * I(3) # cube of $a$ bohrs
# Helium at the center of the box
atoms = [ElementPsp(:He; psp=load_psp("hgh/lda/He-q2"))]
positions = [[1/2, 1/2, 1/2]]
kgrid = [1, 1, 1] # no k-point sampling for an isolated system
Ecut = 30
tol = 1e-8
# dipole moment of a given density (assuming the current geometry)
function dipole(basis, ρ)
rr = [(r[1] - a/2) for r in r_vectors_cart(basis)]
sum(rr .* ρ) * basis.dvol
end;
We first compute the polarizability by finite differences. First compute the dipole moment at rest:
model = model_LDA(lattice, atoms, positions; symmetries=false)
basis = PlaneWaveBasis(model; Ecut, kgrid)
res = self_consistent_field(basis; tol)
μref = dipole(basis, res.ρ)
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -2.770413439866 -0.52 9.0 168ms 2 -2.771692575813 -2.89 -1.32 1.0 89.3ms 3 -2.771714324930 -4.66 -2.44 1.0 89.9ms 4 -2.771714638766 -6.50 -3.14 1.0 93.8ms 5 -2.771714715021 -7.12 -4.27 2.0 106ms 6 -2.771714715144 -9.91 -4.47 1.0 93.6ms 7 -2.771714715248 -9.98 -5.77 1.0 97.6ms 8 -2.771714715250 -11.74 -6.34 2.0 113ms 9 -2.771714715250 + -14.15 -6.36 1.0 108ms 10 -2.771714715250 -14.03 -7.70 1.0 107ms 11 -2.771714715250 + -13.91 -8.41 2.0 121ms
-0.00013457334651682272
Then in a small uniform field:
ε = .01
model_ε = model_LDA(lattice, atoms, positions;
extra_terms=[ExternalFromReal(r -> -ε * (r[1] - a/2))],
symmetries=false)
basis_ε = PlaneWaveBasis(model_ε; Ecut, kgrid)
res_ε = self_consistent_field(basis_ε; tol)
με = dipole(basis_ε, res_ε.ρ)
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -2.770331875483 -0.53 9.0 169ms 2 -2.771767784794 -2.84 -1.30 1.0 88.9ms 3 -2.771801572356 -4.47 -2.64 1.0 98.3ms 4 -2.771802051324 -6.32 -3.98 1.0 95.5ms 5 -2.771802074405 -7.64 -4.68 2.0 114ms 6 -2.771802074467 -10.21 -5.39 1.0 103ms 7 -2.771802074475 -11.09 -5.66 2.0 117ms 8 -2.771802074476 -12.30 -5.96 1.0 103ms 9 -2.771802074476 -12.66 -7.13 1.0 114ms 10 -2.771802074476 + -14.07 -7.24 2.0 130ms 11 -2.771802074476 + -13.83 -7.99 1.0 103ms 12 -2.771802074476 -14.75 -8.53 1.0 109ms
0.017612222131398646
polarizability = (με - μref) / ε
println("Reference dipole: $μref")
println("Displaced dipole: $με")
println("Polarizability : $polarizability")
Reference dipole: -0.00013457334651682272 Displaced dipole: 0.017612222131398646 Polarizability : 1.7746795477915467
The result on more converged grids is very close to published results. For example DOI 10.1039/C8CP03569E quotes 1.65 with LSDA and 1.38 with CCSD(T).
Now we use linear response to compute this analytically; we refer to standard textbooks for the formalism. In the following, $χ_0$ is the independent-particle polarizability, and $K$ the Hartree-exchange-correlation kernel. We denote with $δV_{\rm ext}$ an external perturbing potential (like in this case the uniform electric field). Then: $$ δρ = χ_0 δV = χ_0 (δV_{\rm ext} + K δρ), $$ which implies $$ δρ = (1-χ_0 K)^{-1} χ_0 δV_{\rm ext}. $$ From this we identify the polarizability operator to be $χ = (1-χ_0 K)^{-1} χ_0$. Numerically, we apply $χ$ to $δV = -x$ by solving a linear equation (the Dyson equation) iteratively.
using KrylovKit
# Apply $(1- χ_0 K)$
function dielectric_operator(δρ)
δV = apply_kernel(basis, δρ; res.ρ)
χ0δV = apply_χ0(res, δV)
δρ - χ0δV
end
# `δVext` is the potential from a uniform field interacting with the dielectric dipole
# of the density.
δVext = [-(r[1] - a/2) for r in r_vectors_cart(basis)]
δVext = cat(δVext; dims=4)
# Apply $χ_0$ once to get non-interacting dipole
δρ_nointeract = apply_χ0(res, δVext)
# Solve Dyson equation to get interacting dipole
δρ = linsolve(dielectric_operator, δρ_nointeract, verbosity=3)[1]
println("Non-interacting polarizability: $(dipole(basis, δρ_nointeract))")
println("Interacting polarizability: $(dipole(basis, δρ))")
WARNING: using ChainRulesCore.add!! in module KrylovKit conflicts with an existing identifier. WARNING: using KrylovKit.basis in module ##325 conflicts with an existing identifier. [ Info: GMRES linsolve in iter 1; step 1: normres = 2.493920642759e-01 [ Info: GMRES linsolve in iter 1; step 2: normres = 3.766551036154e-03 [ Info: GMRES linsolve in iter 1; step 3: normres = 2.852762632273e-04 [ Info: GMRES linsolve in iter 1; step 4: normres = 4.694593847233e-06 [ Info: GMRES linsolve in iter 1; step 5: normres = 1.088787643285e-08 [ Info: GMRES linsolve in iter 1; step 6: normres = 6.275052611108e-11 [ Info: GMRES linsolve in iter 1; step 7: normres = 7.187042777764e-13 [ Info: GMRES linsolve in iter 1; finished at step 7: normres = 7.187042777764e-13 [ Info: GMRES linsolve in iter 2; step 1: normres = 1.711137621230e-09 [ Info: GMRES linsolve in iter 2; step 2: normres = 3.356108630396e-11 [ Info: GMRES linsolve in iter 2; step 3: normres = 4.914037638499e-12 [ Info: GMRES linsolve in iter 2; step 4: normres = 6.317995325673e-14 [ Info: GMRES linsolve in iter 2; finished at step 4: normres = 6.317995325673e-14 ┌ Info: GMRES linsolve converged at iteration 2, step 4: │ * norm of residual = 6.295538228293946e-14 └ * number of operations = 13 Non-interacting polarizability: 1.9257125371668444 Interacting polarizability: 1.77365485888704
As expected, the interacting polarizability matches the finite difference result. The non-interacting polarizability is higher.