Polarizability using automatic differentiation
Simple example for computing properties using (forward-mode) automatic differentiation. For a more classical approach and more details about computing polarizabilities, see Polarizability by linear response.
using DFTK
using LinearAlgebra
using ForwardDiff
using PseudoPotentialData
# Construct PlaneWaveBasis given a particular electric field strength
# Again we take the example of a Helium atom.
function make_basis(ε::T; a=10., Ecut=30) where {T}
lattice = T(a) * I(3) # lattice is a cube of ``a`` Bohrs
# Helium at the center of the box
pseudopotentials = PseudoFamily("cp2k.nc.sr.lda.v0_1.semicore.gth")
atoms = [ElementPsp(:He, pseudopotentials)]
positions = [[1/2, 1/2, 1/2]]
model = model_DFT(lattice, atoms, positions;
functionals=[:lda_x, :lda_c_vwn],
extra_terms=[ExternalFromReal(r -> -ε * (r[1] - a/2))],
symmetries=false)
PlaneWaveBasis(model; Ecut, kgrid=[1, 1, 1]) # No k-point sampling on isolated system
end
# dipole moment of a given density (assuming the current geometry)
function dipole(basis, ρ)
@assert isdiag(basis.model.lattice)
a = basis.model.lattice[1, 1]
rr = [a * (r[1] - 1/2) for r in r_vectors(basis)]
sum(rr .* ρ) * basis.dvol
end
# Function to compute the dipole for a given field strength
function compute_dipole(ε; tol=1e-8, kwargs...)
scfres = self_consistent_field(make_basis(ε; kwargs...); tol)
dipole(scfres.basis, scfres.ρ)
end;
With this in place we can compute the polarizability from finite differences (just like in the previous example):
polarizability_fd = let
ε = 0.01
(compute_dipole(ε) - compute_dipole(0.0)) / ε
end
1.773558041225318
We do the same thing using automatic differentiation. Under the hood this uses custom rules to implicitly differentiate through the self-consistent field fixed-point problem. This leads to a density-functional perturbation theory problem, which is automatically set up and solved in the background.
polarizability = ForwardDiff.derivative(compute_dipole, 0.0)
println()
println("Polarizability via ForwardDiff: $polarizability")
println("Polarizability via finite difference: $polarizability_fd")
n Energy log10(ΔE) log10(Δρ) Diag Δtime
--- --------------- --------- --------- ---- ------
1 -2.770715706660 -0.52 8.0 154ms
2 -2.772053532041 -2.87 -1.31 1.0 130ms
3 -2.772083064229 -4.53 -2.53 1.0 117ms
4 -2.772083368468 -6.52 -3.40 1.0 145ms
5 -2.772083415721 -7.33 -3.88 2.0 120ms
6 -2.772083417765 -8.69 -5.27 1.0 109ms
7 -2.772083417806 -10.39 -5.21 2.0 123ms
8 -2.772083417810 -11.39 -5.60 1.0 114ms
9 -2.772083417811 -12.34 -6.98 1.0 113ms
10 -2.772083417811 -13.34 -6.88 2.0 125ms
11 -2.772083417811 -13.88 -7.44 1.0 115ms
12 -2.772083417811 + -14.07 -7.67 1.0 124ms
13 -2.772083417811 -14.18 -8.00 1.0 154ms
Solving response problem
[ Info: GMRES linsolve starts with norm of residual = 4.19e+00
[ Info: GMRES linsolve in iteration 1; step 1: normres = 2.49e-01
[ Info: GMRES linsolve in iteration 1; step 2: normres = 3.76e-03
[ Info: GMRES linsolve in iteration 1; step 3: normres = 2.84e-04
[ Info: GMRES linsolve in iteration 1; step 4: normres = 4.67e-06
[ Info: GMRES linsolve in iteration 1; step 5: normres = 1.08e-08
[ Info: GMRES linsolve in iteration 1; step 6: normres = 2.35e-08
┌ Info: GMRES linsolve converged at iteration 2, step 1:
│ * norm of residual = 1.79e-09
└ * number of operations = 9
Polarizability via ForwardDiff: 1.7725349973039477
Polarizability via finite difference: 1.773558041225318