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

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.770801903669                   -0.52    9.0    167ms
  2   -2.772059775620       -2.90       -1.32    1.0    189ms
  3   -2.772083091457       -4.63       -2.47    1.0   91.5ms
  4   -2.772083342248       -6.60       -3.18    1.0    124ms
  5   -2.772083414562       -7.14       -3.77    2.0    110ms
  6   -2.772083417700       -8.50       -4.79    1.0    100ms
  7   -2.772083417809       -9.96       -5.35    2.0    118ms
  8   -2.772083417810      -12.12       -5.85    1.0    105ms
  9   -2.772083417811      -12.29       -7.16    2.0    123ms
 10   -2.772083417811   +  -14.27       -7.51    1.0    105ms
 11   -2.772083417811      -14.21       -7.93    1.0    108ms
 12   -2.772083417811      -14.51       -9.68    1.0    106ms
Solving response problem
Iter  Restart  Krydim  log10(res)  avg(CG)  Δtime   Comment
----  -------  ------  ----------  -------  ------  ---------------
                                      15.0  70.9ms  Non-interacting
   1        0       1       -0.60     12.0   244ms  
   2        0       2       -2.42     10.0   122ms  
   3        0       3       -3.55      7.0   142ms  
   4        0       4       -5.33      6.0   121ms  
   5        0       5       -7.97      4.0   101ms  
   6        0       6      -10.26      1.0  90.2ms  
   7        1       1       -9.21     14.0   232ms  Restart
   8        1       2      -10.79      1.0  93.1ms  
                                      14.0   133ms  Final orbitals

Polarizability via ForwardDiff:       1.772534970624385
Polarizability via finite difference: 1.7735582266998096