{ "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\n", "--- --------------- --------- --------- ----\n", " 1 -2.770831422053 -0.53 8.0\n", " 2 -2.772140124772 -2.88 -1.31 1.0\n", " 3 -2.772170111817 -4.52 -2.58 1.0\n", " 4 -2.772170715344 -6.22 -3.71 2.0\n", " 5 -2.772170722345 -8.15 -4.11 2.0\n", " 6 -2.772170723002 -9.18 -4.99 1.0\n", " 7 -2.772170723014 -10.93 -5.77 1.0\n", " 8 -2.772170723015 -11.94 -6.30 2.0\n", " 9 -2.772170723015 -14.03 -6.69 1.0\n", " 10 -2.772170723015 -14.45 -7.09 1.0\n", " 11 -2.772170723015 + -15.35 -7.90 2.0\n", " 12 -2.772170723015 -13.92 -8.55 1.0\n", "n Energy log10(ΔE) log10(Δρ) Diag\n", "--- --------------- --------- --------- ----\n", " 1 -2.770789244638 -0.52 9.0\n", " 2 -2.772060985538 -2.90 -1.32 1.0\n", " 3 -2.772083109637 -4.66 -2.47 1.0\n", " 4 -2.772083345932 -6.63 -3.17 1.0\n", " 5 -2.772083416701 -7.15 -3.98 2.0\n", " 6 -2.772083417755 -8.98 -4.99 1.0\n", " 7 -2.772083417811 -10.25 -5.66 2.0\n", " 8 -2.772083417811 -13.94 -6.21 1.0\n", " 9 -2.772083417811 -13.16 -6.94 2.0\n", " 10 -2.772083417811 + -13.97 -6.85 1.0\n", " 11 -2.772083417811 -13.92 -7.41 1.0\n", " 12 -2.772083417811 + -14.35 -7.54 1.0\n", " 13 -2.772083417811 + -14.51 -7.80 1.0\n", " 14 -2.772083417811 -14.88 -7.92 1.0\n", " 15 -2.772083417811 -14.07 -8.47 1.0\n" ] }, { "output_type": "execute_result", "data": { "text/plain": "1.773558037674741" }, "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\n", "--- --------------- --------- --------- ----\n", " 1 -2.770771102104 -0.52 9.0\n", " 2 -2.772060864673 -2.89 -1.32 1.0\n", " 3 -2.772083027561 -4.65 -2.44 1.0\n", " 4 -2.772083338213 -6.51 -3.13 1.0\n", " 5 -2.772083417490 -7.10 -4.22 2.0\n", " 6 -2.772083417702 -9.68 -4.55 1.0\n", " 7 -2.772083417809 -9.97 -5.69 1.0\n", " 8 -2.772083417811 -11.93 -6.07 2.0\n", " 9 -2.772083417811 -13.42 -6.43 1.0\n", " 10 -2.772083417811 -13.89 -8.25 1.0\n", "\n", "Polarizability via ForwardDiff: 1.7725349319705532\n", "Polarizability via finite difference: 1.773558037674741\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.8.4" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.4", "language": "julia" } }, "nbformat": 4 }