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.7735579636506626

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.770663446679                   -0.53    8.0    217ms
  2   -2.772053224810       -2.86       -1.30    1.0   92.3ms
  3   -2.772082722982       -4.53       -2.63    1.0    118ms
  4   -2.772083413436       -6.16       -3.89    2.0    114ms
  5   -2.772083417374       -8.40       -4.23    2.0    116ms
  6   -2.772083417798       -9.37       -5.32    1.0    102ms
  7   -2.772083417810      -10.89       -5.90    2.0    122ms
  8   -2.772083417811      -12.77       -6.56    1.0    104ms
  9   -2.772083417811      -14.65       -7.30    2.0    123ms
 10   -2.772083417811      -14.31       -7.81    1.0    189ms
 11   -2.772083417811   +  -14.24       -8.82    2.0    871ms
Solving response problem
Iter  Restart  Krydim  log10(res)  avg(CG)  Δtime   Comment
----  -------  ------  ----------  -------  ------  ---------------
                                      13.0  65.5ms  Non-interacting
   1        0       1       -0.60     11.0   207ms  
   2        0       2       -2.42      8.0   146ms  
   3        0       3       -3.55      6.0   124ms  
   4        0       4       -5.33      5.0   108ms  
   5        0       5       -7.93      2.0  95.6ms  
   6        0       6      -10.45      1.0  91.4ms  
   7        1       1       -8.13     14.0   241ms  Restart
   8        1       2      -10.41      1.0  90.7ms  
                                      13.0   138ms  Final orbitals

Polarizability via ForwardDiff:       1.772534970641576
Polarizability via finite difference: 1.7735579636506626