{ "cells": [ { "cell_type": "markdown", "source": [ "# Polarizability using automatic differentiation\n", "\n", "Simple example for computing properties using (forward-mode)\n", "automatic differentiation.\n", "For a more classical approach and more details about computing polarizabilities,\n", "see Polarizability by linear response." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using LinearAlgebra\n", "using ForwardDiff\n", "\n", "# Construct PlaneWaveBasis given a particular electric field strength\n", "# Again we take the example of a Helium atom.\n", "function make_basis(ε::T; a=10., Ecut=30) where {T}\n", " lattice=T(a) * I(3) # lattice is a cube of $a$ Bohrs\n", " # Helium at the center of the box\n", " atoms = [ElementPsp(:He, psp=load_psp(\"hgh/lda/He-q2\"))]\n", " positions = [[1/2, 1/2, 1/2]]\n", "\n", " model = model_DFT(lattice, atoms, positions, [:lda_x, :lda_c_vwn];\n", " extra_terms=[ExternalFromReal(r -> -ε * (r[1] - a/2))],\n", " symmetries=false)\n", " PlaneWaveBasis(model; Ecut, kgrid=[1, 1, 1]) # No k-point sampling on isolated system\n", "end\n", "\n", "# dipole moment of a given density (assuming the current geometry)\n", "function dipole(basis, ρ)\n", " @assert isdiag(basis.model.lattice)\n", " a = basis.model.lattice[1, 1]\n", " rr = [a * (r[1] - 1/2) for r in r_vectors(basis)]\n", " sum(rr .* ρ) * basis.dvol\n", "end\n", "\n", "# Function to compute the dipole for a given field strength\n", "function compute_dipole(ε; tol=1e-8, kwargs...)\n", " scfres = self_consistent_field(make_basis(ε; kwargs...); tol)\n", " dipole(scfres.basis, scfres.ρ)\n", "end;" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "With this in place we can compute the polarizability from finite differences\n", "(just like in the previous example):" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) Diag Δtime\n", "--- --------------- --------- --------- ---- ------\n", " 1 -2.770882126745 -0.52 9.0 \n", " 2 -2.772147984284 -2.90 -1.32 1.0 158ms\n", " 3 -2.772170271702 -4.65 -2.46 1.0 138ms\n", " 4 -2.772170653526 -6.42 -3.17 1.0 158ms\n", " 5 -2.772170722681 -7.16 -4.14 2.0 159ms\n", " 6 -2.772170722846 -9.78 -4.31 1.0 143ms\n", " 7 -2.772170723011 -9.78 -5.32 1.0 156ms\n", " 8 -2.772170723015 -11.42 -5.95 2.0 165ms\n", " 9 -2.772170723015 -13.03 -6.57 1.0 160ms\n", " 10 -2.772170723015 + -14.21 -6.50 2.0 176ms\n", " 11 -2.772170723015 -13.97 -7.59 1.0 178ms\n", " 12 -2.772170723015 + -14.75 -7.82 2.0 167ms\n", " 13 -2.772170723015 + -14.65 -8.42 1.0 171ms\n", "n Energy log10(ΔE) log10(Δρ) Diag Δtime\n", "--- --------------- --------- --------- ---- ------\n", " 1 -2.770804760980 -0.52 9.0 \n", " 2 -2.772061412731 -2.90 -1.32 1.0 157ms\n", " 3 -2.772083168551 -4.66 -2.48 1.0 138ms\n", " 4 -2.772083351796 -6.74 -3.19 1.0 153ms\n", " 5 -2.772083415788 -7.19 -3.87 2.0 159ms\n", " 6 -2.772083417695 -8.72 -4.70 1.0 142ms\n", " 7 -2.772083417796 -10.00 -4.94 2.0 170ms\n", " 8 -2.772083417810 -10.85 -5.73 1.0 148ms\n", " 9 -2.772083417811 -12.23 -6.15 2.0 187ms\n", " 10 -2.772083417811 -13.17 -7.64 1.0 150ms\n", " 11 -2.772083417811 -14.57 -8.03 2.0 186ms\n" ] }, { "output_type": "execute_result", "data": { "text/plain": "1.7735579529824077" }, "metadata": {}, "execution_count": 2 } ], "cell_type": "code", "source": [ "polarizability_fd = let\n", " ε = 0.01\n", " (compute_dipole(ε) - compute_dipole(0.0)) / ε\n", "end" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "We do the same thing using automatic differentiation. Under the hood this uses\n", "custom rules to implicitly differentiate through the self-consistent\n", "field fixed-point problem." ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) Diag Δtime\n", "--- --------------- --------- --------- ---- ------\n", " 1 -2.770738982964 -0.52 8.0 \n", " 2 -2.772057858034 -2.88 -1.32 1.0 136ms\n", " 3 -2.772083163774 -4.60 -2.48 1.0 137ms\n", " 4 -2.772083350523 -6.73 -3.21 1.0 173ms\n", " 5 -2.772083414090 -7.20 -3.76 2.0 158ms\n", " 6 -2.772083417661 -8.45 -4.69 1.0 157ms\n", " 7 -2.772083417806 -9.84 -5.12 2.0 161ms\n", " 8 -2.772083417810 -11.33 -6.07 1.0 159ms\n", " 9 -2.772083417811 -12.61 -6.87 2.0 174ms\n", " 10 -2.772083417811 -14.10 -7.69 1.0 161ms\n", " 11 -2.772083417811 -13.94 -8.15 2.0 168ms\n", "\n", "Polarizability via ForwardDiff: 1.772534962791999\n", "Polarizability via finite difference: 1.7735579529824077\n" ] } ], "cell_type": "code", "source": [ "polarizability = ForwardDiff.derivative(compute_dipole, 0.0)\n", "println()\n", "println(\"Polarizability via ForwardDiff: $polarizability\")\n", "println(\"Polarizability via finite difference: $polarizability_fd\")" ], "metadata": {}, "execution_count": 3 } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.9.0" }, "kernelspec": { "name": "julia-1.9", "display_name": "Julia 1.9.0", "language": "julia" } }, "nbformat": 4 }