Polarizability by linear response

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,

In [1]:
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, ρ)
    rr = [a * (r[1] - 1/2) for r in r_vectors(basis)]
    d = sum(rr .* ρ) * basis.dvol
end;

Polarizability by finite differences

We first compute the polarizability by finite differences. First compute the dipole moment at rest:

In [2]:
model = model_LDA(lattice, atoms; symmetries=false)
basis = PlaneWaveBasis(model, Ecut; kgrid=kgrid)
res = self_consistent_field(basis, tol=tol)
μref = dipole(basis, res.ρ)
n     Energy            Eₙ-Eₙ₋₁     ρout-ρin   Diag
---   ---------------   ---------   --------   ----
  1   -2.769904414614         NaN   2.99e-01    9.0 
  2   -2.771189478564   -1.29e-03   4.74e-02    1.0 
  3   -2.771212028375   -2.25e-05   4.51e-03    2.0 
  4   -2.771212981445   -9.53e-07   1.10e-04    2.0 
  5   -2.771212981733   -2.88e-10   1.51e-05    2.0 
Out[2]:
-0.00013510287592493662

Then in a small uniform field:

In [3]:
ε = .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(basis_ε, res_ε.ρ)
n     Energy            Eₙ-Eₙ₋₁     ρout-ρin   Diag
---   ---------------   ---------   --------   ----
  1   -2.769991270361         NaN   2.99e-01    8.0 
  2   -2.771275671868   -1.28e-03   4.77e-02    1.0 
  3   -2.771299717229   -2.40e-05   4.11e-03    2.0 
  4   -2.771300361667   -6.44e-07   6.98e-05    2.0 
  5   -2.771300362343   -6.76e-10   1.40e-05    2.0 
Out[3]:
0.01761886874200919
In [4]:
polarizability = (με - μref) / ε

println("Reference dipole:  $μref")
println("Displaced dipole:  $με")
println("Polarizability :   $polarizability")
Reference dipole:  -0.00013510287592493662
Displaced dipole:  0.01761886874200919
Polarizability :   1.7753971617934128

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: $$ d\rho = \chi_0 dV = \chi_0 (dV_{\rm ext} + K d\rho), $$ which implies $$ d\rho = (1-\chi_0 K)^-1 \chi_0 dV_{\rm ext}. $$ 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.

In [5]:
using KrylovKit

# Apply (1- χ0 K)
function dielectric_operator()
    dv = apply_kernel(basis, ; ρ=res.ρ)
    χ0dv = apply_χ0(res.ham, res.ψ, res.εF, res.eigenvalues, dv)
     - χ0dv
end

# dVext is the potential from a uniform field interacting with the dielectric dipole
# of the density.
dVext = [-a * (r[1] - 1/2) for r in r_vectors(basis)]
dVext = cat(dVext; dims=4)

# Apply χ0 once to get non-interacting dipole
dρ_nointeract = apply_χ0(res.ham, res.ψ, res.εF, res.eigenvalues, dVext)

# Solve Dyson equation to get interacting dipole
 = linsolve(dielectric_operator, dρ_nointeract, verbosity=3)[1]

println("Non-interacting polarizability: $(dipole(basis, dρ_nointeract))")
println("Interacting polarizability:     $(dipole(basis, ))")
WARNING: using KrylovKit.basis in module ##278 conflicts with an existing identifier.
┌ Info: GMRES linsolve in iter 1; step 1: normres = 2.490073984786e-01
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:55
┌ Info: GMRES linsolve in iter 1; step 2: normres = 4.823172401830e-03
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 1; step 3: normres = 8.978070085058e-04
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 1; step 4: normres = 7.660510694749e-06
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 1; step 5: normres = 6.898680811602e-07
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 1; step 6: normres = 1.944534643258e-09
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 1; step 7: normres = 3.359263353795e-11
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 1; step 8: normres = 3.025629400497e-12
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 1; finished at step 8: normres = 3.025629400497e-12
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:96
┌ Info: GMRES linsolve in iter 2; step 1: normres = 1.402304215181e-06
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:55
┌ Info: GMRES linsolve in iter 2; step 2: normres = 3.307084868319e-08
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 2; step 3: normres = 1.342885281636e-09
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 2; step 4: normres = 2.511365741322e-11
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 2; step 5: normres = 1.723203596002e-13
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 2; finished at step 5: normres = 1.723203596002e-13
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:96
┌ Info: GMRES linsolve in iter 3; step 1: normres = 1.880301919106e-11
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:55
┌ Info: GMRES linsolve in iter 3; step 2: normres = 1.629217282114e-12
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:89
┌ Info: GMRES linsolve in iter 3; finished at step 2: normres = 1.629217282114e-12
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:96
┌ Info: GMRES linsolve converged at iteration 3, step 2:
│ *  norm of residual = 1.6178211483524532e-12
│ *  number of operations = 17
└ @ KrylovKit /home/runner/.julia/packages/KrylovKit/YPiz7/src/linsolve/gmres.jl:127
Non-interacting polarizability: 1.9300222967209641
Interacting polarizability:     1.7765407914986935

As expected, the interacting polarizability matches the finite difference result. The non-interacting polarizability is higher.