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.770383275175 -0.52 9.0 2 -2.771690376929 -2.88 -1.33 1.0 184ms 3 -2.771714227464 -4.62 -2.42 1.0 195ms 4 -2.771714621352 -6.40 -3.10 1.0 236ms 5 -2.771714714971 -7.03 -4.25 2.0 207ms 6 -2.771714715146 -9.76 -4.56 1.0 191ms 7 -2.771714715248 -9.99 -5.90 1.0 228ms 8 -2.771714715250 -11.92 -6.35 2.0 236ms 9 -2.771714715250 -13.51 -6.50 1.0 234ms 10 -2.771714715250 -13.40 -7.75 1.0 236ms 11 -2.771714715250 -14.88 -8.53 2.0 277ms
-0.00013457367159238574
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.770471813681 -0.53 8.0 2 -2.771769776472 -2.89 -1.31 1.0 238ms 3 -2.771801185881 -4.50 -2.56 1.0 187ms 4 -2.771802065966 -6.06 -3.70 2.0 218ms 5 -2.771802073871 -8.10 -4.12 2.0 250ms 6 -2.771802074459 -9.23 -5.32 1.0 189ms 7 -2.771802074476 -10.78 -6.14 2.0 254ms 8 -2.771802074476 -13.30 -6.45 1.0 210ms 9 -2.771802074476 -13.37 -7.18 2.0 248ms 10 -2.771802074476 + -14.10 -7.63 1.0 243ms 11 -2.771802074476 -15.05 -8.22 2.0 251ms
0.017612221419371757
polarizability = (με - μref) / ε
println("Reference dipole: $μref")
println("Displaced dipole: $με")
println("Polarizability : $polarizability")
Reference dipole: -0.00013457367159238574 Displaced dipole: 0.017612221419371757 Polarizability : 1.774679509096414
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 KrylovKit.basis in module ##352 conflicts with an existing identifier. [ Info: GMRES linsolve in iter 1; step 1: normres = 2.493920806915e-01 [ Info: GMRES linsolve in iter 1; step 2: normres = 3.766551546892e-03 [ Info: GMRES linsolve in iter 1; step 3: normres = 2.852768076049e-04 [ Info: GMRES linsolve in iter 1; step 4: normres = 4.694603439891e-06 [ Info: GMRES linsolve in iter 1; step 5: normres = 1.088788386421e-08 [ Info: GMRES linsolve in iter 1; step 6: normres = 6.289880553836e-11 [ Info: GMRES linsolve in iter 1; step 7: normres = 1.171550244880e-12 [ Info: GMRES linsolve in iter 1; finished at step 7: normres = 1.171550244880e-12 [ Info: GMRES linsolve in iter 2; step 1: normres = 9.456243833546e-11 [ Info: GMRES linsolve in iter 2; step 2: normres = 1.500959869909e-11 [ Info: GMRES linsolve in iter 2; step 3: normres = 3.169443494038e-13 [ Info: GMRES linsolve in iter 2; finished at step 3: normres = 3.169443494038e-13 ┌ Info: GMRES linsolve converged at iteration 2, step 3: │ * norm of residual = 3.1690805851827247e-13 └ * number of operations = 12 Non-interacting polarizability: 1.9257125406109241 Interacting polarizability: 1.7736548665321363
As expected, the interacting polarizability matches the finite difference result. The non-interacting polarizability is higher.