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

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.770737574769                   -0.53    9.0    171ms
  2   -2.772051715666       -2.88       -1.31    1.0   99.0ms
  3   -2.772082838812       -4.51       -2.59    1.0    200ms
  4   -2.772083379678       -6.27       -3.93    1.0    102ms
  5   -2.772083417249       -7.43       -4.17    2.0    154ms
  6   -2.772083417775       -9.28       -5.14    1.0    106ms
  7   -2.772083417808      -10.48       -5.37    2.0    130ms
  8   -2.772083417810      -11.67       -5.82    1.0    111ms
  9   -2.772083417811      -12.33       -6.54    2.0    131ms
 10   -2.772083417811      -14.10       -7.22    1.0    115ms
 11   -2.772083417811      -14.35       -8.00    1.0    116ms
Solving response problem
Iter  Restart  Krydim  log10(res)  avg(CG)  Δtime   Comment
----  -------  ------  ----------  -------  ------  ---------------
                                      13.0  59.6ms  Non-interacting
   1        0       1       -0.60     10.0   178ms
   2        0       2       -2.42      8.0   131ms
   3        0       3       -3.55      5.0   103ms
   4        0       4       -5.33      4.0  98.8ms
   5        0       5       -7.40      1.0  87.7ms
   6        0       6      -10.28      1.0  87.8ms
   7        1       1       -7.42     12.0   221ms  Restart
   8        1       2       -9.04      1.0  87.5ms
                                      13.0   122ms  Final orbitals

Polarizability via ForwardDiff:       1.7725349928572163
Polarizability via finite difference: 1.7735579478260433