Polarizability by linear response
We compute the polarizability of a Helium atom. The polarizability is defined as the change in dipole moment
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
He = ElementPsp(:He, psp=load_psp("hgh/lda/He-q2"))
atoms = [He => [[1/2; 1/2; 1/2]]] # Helium at the center of the box
kgrid = [1, 1, 1] # no kpoint sampling for an isolated system
Ecut = 30
tol = 1e-8
# dipole moment of a given density (assuming the current geometry)
function dipole(ρ)
basis = ρ.basis
rr = [a * (r[1] - 1/2) for r in r_vectors(basis)]
d = sum(rr .* ρ.real) * basis.model.unit_cell_volume / prod(basis.fft_size)
end;
Polarizability by finite differences
We first compute the polarizability by finite differences. First compute the dipole moment at rest:
model = model_LDA(lattice, atoms; symmetries=false)
basis = PlaneWaveBasis(model, Ecut; kgrid=kgrid)
res = self_consistent_field(basis, tol=tol)
μref = dipole(res.ρ)
-0.0001347981222375538
Then in a small uniform field:
ε = .01
model_ε = model_LDA(lattice, atoms; extra_terms=[ExternalFromReal(r -> -ε * (r[1] - a/2))],
symmetries=false)
basis_ε = PlaneWaveBasis(model_ε, Ecut; kgrid=kgrid)
res_ε = self_consistent_field(basis_ε, tol=tol)
με = dipole(res_ε.ρ)
0.017614932944754965
polarizability = (με - μref) / ε
println("Reference dipole: $μref")
println("Displaced dipole: $με")
println("Polarizability : $polarizability")
Reference dipole: -0.0001347981222375538 Displaced dipole: 0.017614932944754965 Polarizability : 1.774973106699252
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).
Polarizability by linear response
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 $dV_{\rm ext}$ an external perturbing potential (like in this case the uniform electric field). Then:
which implies
From this we identify the polarizability operator to be $\chi = (1-\chi_0 K)^-1 \chi_0$. Numerically, we apply $\chi$ to $dV = -x$ by solving a linear equation (the Dyson equation) iteratively.
using KrylovKit
# KrylovKit cannot deal with the density as a 3D array, so we need `vec` to vectorize
# and `devec` to "unvectorize"
devec(arr) = from_real(basis, reshape(arr, size(res.ρ.real)))
# Apply (1- χ0 K) to a vectorized dρ
function dielectric_operator(dρ)
dρ = devec(dρ)
dv = apply_kernel(basis, dρ; ρ=res.ρ)
χ0dv = apply_χ0(res.ham, res.ψ, res.εF, res.eigenvalues, dv...)[1]
vec((dρ - χ0dv).real)
end
# dVext is the potential from a uniform field interacting with the dielectric dipole
# of the density.
dVext = from_real(basis, [-a * (r[1] - 1/2) for r in r_vectors(basis)])
# Apply χ0 once to get non-interacting dipole
dρ_nointeract = apply_χ0(res.ham, res.ψ, res.εF, res.eigenvalues, dVext)[1]
# Solve Dyson equation to get interacting dipole
dρ = devec(linsolve(dielectric_operator, vec(dρ_nointeract.real), verbosity=3)[1])
println("Non-interacting polarizability: $(dipole(dρ_nointeract))")
println("Interacting polarizability: $(dipole(dρ))")
WARNING: using KrylovKit.basis in module ex-polarizability conflicts with an existing identifier. [ Info: GMRES linsolve in iter 1; step 1: normres = 2.492058630690e-01 [ Info: GMRES linsolve in iter 1; step 2: normres = 3.761903585717e-03 [ Info: GMRES linsolve in iter 1; step 3: normres = 2.921940883384e-04 [ Info: GMRES linsolve in iter 1; step 4: normres = 4.692162513740e-06 [ Info: GMRES linsolve in iter 1; step 5: normres = 7.064413607669e-08 [ Info: GMRES linsolve in iter 1; step 6: normres = 1.838445710054e-09 [ Info: GMRES linsolve in iter 1; step 7: normres = 1.089999609407e-11 [ Info: GMRES linsolve in iter 1; step 8: normres = 1.094629519445e-13 [ Info: GMRES linsolve in iter 1; finised at step 8: normres = 1.094629519445e-13 [ Info: GMRES linsolve in iter 2; step 1: normres = 7.917881146464e-07 [ Info: GMRES linsolve in iter 2; step 2: normres = 9.433743030233e-08 [ Info: GMRES linsolve in iter 2; step 3: normres = 1.729892120609e-09 [ Info: GMRES linsolve in iter 2; step 4: normres = 2.849145695613e-11 [ Info: GMRES linsolve in iter 2; step 5: normres = 2.076698420555e-13 [ Info: GMRES linsolve in iter 2; finised at step 5: normres = 2.076698420555e-13 [ Info: GMRES linsolve in iter 3; step 1: normres = 4.030951490743e-11 [ Info: GMRES linsolve in iter 3; step 2: normres = 5.092225758383e-12 [ Info: GMRES linsolve in iter 3; step 3: normres = 4.556546842314e-14 [ Info: GMRES linsolve in iter 3; finised at step 3: normres = 4.556546842314e-14 ┌ Info: GMRES linsolve converged at iteration 3, step 3: │ * norm of residual = 4.603744423650452e-14 └ * number of operations = 18 Non-interacting polarizability: 1.9263134288623363 Interacting polarizability: 1.7740881447430417
As expected, the interacting polarizability matches the finite difference result. The non-interacting polarizability is higher.