We compute the polarizability of a Helium atom. The polarizability is defined as the change in dipole moment $$ \mu = \int 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=tol)
μref = dipole(basis, res.ρ)
n Energy log10(ΔE) log10(Δρ) Diag --- --------------- --------- --------- ---- 1 -2.770368237156 -0.53 8.0 2 -2.771674719268 -2.88 -1.29 1.0 3 -2.771714591553 -4.40 -2.84 2.0 4 -2.771714713326 -6.91 -3.69 2.0 5 -2.771714715212 -8.72 -4.45 2.0
-0.00013425080761394765
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=tol)
με = dipole(basis_ε, res_ε.ρ)
n Energy log10(ΔE) log10(Δρ) Diag --- --------------- --------- --------- ---- 1 -2.770487607571 -0.52 9.0 2 -2.771778814203 -2.89 -1.32 1.0 3 -2.771801233539 -4.65 -2.36 2.0 4 -2.771802073966 -6.08 -4.04 2.0 5 -2.771802074463 -9.30 -4.98 2.0
0.017613515075492962
polarizability = (με - μref) / ε
println("Reference dipole: $μref")
println("Displaced dipole: $με")
println("Polarizability : $polarizability")
Reference dipole: -0.00013425080761394765 Displaced dipole: 0.017613515075492962 Polarizability : 1.774776588310691
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, $\chi_0$ is the independent-particle polarizability, and $K$ the Hartree-exchange-correlation kernel. We denote with $\delta V_{\rm ext}$ an external perturbing potential (like in this case the uniform electric field). Then: $$ \delta\rho = \chi_0 \delta V = \chi_0 (\delta V_{\rm ext} + K \delta\rho), $$ which implies $$ \delta\rho = (1-\chi_0 K)^{-1} \chi_0 \delta V_{\rm ext}. $$ From this we identify the polarizability operator to be $\chi = (1-\chi_0 K)^{-1} \chi_0$. Numerically, we apply $\chi$ to $\delta 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.493960994241e-01 [ Info: GMRES linsolve in iter 1; step 2: normres = 3.766366718930e-03 [ Info: GMRES linsolve in iter 1; step 3: normres = 2.845722327751e-04 [ Info: GMRES linsolve in iter 1; step 4: normres = 4.694835823521e-06 [ Info: GMRES linsolve in iter 1; step 5: normres = 1.093071004865e-08 [ Info: GMRES linsolve in iter 1; step 6: normres = 2.074071454657e-10 [ Info: GMRES linsolve in iter 1; step 7: normres = 1.075635404455e-11 [ Info: GMRES linsolve in iter 1; step 8: normres = 1.123910552290e-13 [ Info: GMRES linsolve in iter 1; finished at step 8: normres = 1.123910552290e-13 [ Info: GMRES linsolve in iter 2; step 1: normres = 8.973332962649e-11 [ Info: GMRES linsolve in iter 2; step 2: normres = 7.799012026293e-12 [ Info: GMRES linsolve in iter 2; step 3: normres = 3.533387081846e-13 [ Info: GMRES linsolve in iter 2; finished at step 3: normres = 3.533387081846e-13 ┌ Info: GMRES linsolve converged at iteration 2, step 3: │ * norm of residual = 3.5335749639422155e-13 └ * number of operations = 13 Non-interacting polarizability: 1.9257641897280278 Interacting polarizability: 1.7737031542947217
As expected, the interacting polarizability matches the finite difference result. The non-interacting polarizability is higher.