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

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.770762203769                   -0.52    9.0    227ms
  2   -2.772057547325       -2.89       -1.32    1.0   96.7ms
  3   -2.772083124392       -4.59       -2.49    1.0    123ms
  4   -2.772083345141       -6.66       -3.22    1.0    103ms
  5   -2.772083413973       -7.16       -3.73    2.0    120ms
  6   -2.772083417621       -8.44       -4.60    1.0    106ms
  7   -2.772083417806       -9.73       -5.13    2.0    125ms
  8   -2.772083417810      -11.42       -5.66    1.0    115ms
  9   -2.772083417811      -12.06       -6.13    2.0    133ms
 10   -2.772083417811      -13.16       -7.46    1.0    123ms
 11   -2.772083417811      -14.65       -7.66    2.0    128ms
 12   -2.772083417811   +  -14.45       -8.15    1.0    114ms
Solving response problem
Iter  Restart  Krydim  log10(res)  avg(CG)  Δtime   Comment
----  -------  ------  ----------  -------  ------  ---------------
                                      13.0   215ms  Non-interacting
   1        0       1       -0.60     10.0   634ms
   2        0       2       -2.42      8.0   118ms
   3        0       3       -3.55      5.0   102ms
   4        0       4       -5.33      5.0   105ms
   5        0       5       -7.53      1.0  84.8ms
   6        0       6      -10.09      1.0  87.9ms
   7        1       1       -7.31     13.0   219ms  Restart
   8        1       2       -8.89      1.0  86.6ms
                                      13.0   127ms  Final orbitals

Polarizability via ForwardDiff:       1.7725349682180154
Polarizability via finite difference: 1.77355795757127