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
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -2.770859815902 -0.53 9.0 177ms 2 -2.772140634348 -2.89 -1.31 1.0 114ms 3 -2.772169899236 -4.53 -2.57 1.0 110ms 4 -2.772170706547 -6.09 -3.53 2.0 127ms 5 -2.772170722269 -7.80 -4.12 1.0 114ms 6 -2.772170722951 -9.17 -4.93 1.0 161ms 7 -2.772170723014 -10.20 -5.41 2.0 496ms 8 -2.772170723015 -12.04 -6.25 1.0 114ms 9 -2.772170723015 -14.31 -6.79 1.0 139ms 10 -2.772170723015 -14.12 -7.47 1.0 115ms 11 -2.772170723015 -14.07 -8.22 2.0 136ms n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -2.770755627493 -0.53 8.0 155ms 2 -2.772053846938 -2.89 -1.31 1.0 109ms 3 -2.772082926712 -4.54 -2.55 1.0 108ms 4 -2.772083368510 -6.35 -3.45 1.0 111ms 5 -2.772083416412 -7.32 -3.96 2.0 127ms 6 -2.772083417755 -8.87 -5.19 1.0 117ms 7 -2.772083417807 -10.28 -5.28 2.0 131ms 8 -2.772083417811 -11.50 -6.30 1.0 121ms 9 -2.772083417811 -13.23 -6.56 2.0 539ms 10 -2.772083417811 -14.18 -7.40 1.0 114ms 11 -2.772083417811 + -14.10 -8.00 1.0 116ms 12 -2.772083417811 -13.82 -8.66 2.0 154ms
1.7735579730013822
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.770657324401 -0.53 9.0 180ms 2 -2.772053993343 -2.85 -1.30 1.0 106ms 3 -2.772082856401 -4.54 -2.65 1.0 110ms 4 -2.772083415438 -6.25 -4.03 2.0 127ms 5 -2.772083417616 -8.66 -4.45 2.0 128ms 6 -2.772083417803 -9.73 -5.42 1.0 600ms 7 -2.772083417810 -11.12 -5.99 2.0 148ms 8 -2.772083417811 -13.01 -6.50 1.0 152ms 9 -2.772083417811 -13.64 -7.01 2.0 145ms 10 -2.772083417811 + -15.35 -7.92 1.0 116ms 11 -2.772083417811 + -14.35 -8.56 2.0 130ms 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 converged at iteration 1, step 6: │ * norm of residual = 7.68e-10 └ * number of operations = 8 Polarizability via ForwardDiff: 1.7725349649275122 Polarizability via finite difference: 1.7735579730013822