{ "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.770805509700 -0.52 9.0 \n", " 2 -2.772145785646 -2.87 -1.32 1.0 163ms\n", " 3 -2.772170260912 -4.61 -2.46 1.0 166ms\n", " 4 -2.772170638205 -6.42 -3.15 1.0 214ms\n", " 5 -2.772170722585 -7.07 -4.08 2.0 193ms\n", " 6 -2.772170722745 -9.80 -4.21 1.0 172ms\n", " 7 -2.772170722984 -9.62 -4.70 1.0 203ms\n", " 8 -2.772170723009 -10.60 -5.09 1.0 178ms\n", " 9 -2.772170723015 -11.24 -6.22 1.0 200ms\n", " 10 -2.772170723015 -12.82 -6.19 2.0 219ms\n", " 11 -2.772170723015 -13.24 -6.84 1.0 183ms\n", " 12 -2.772170723015 -14.01 -7.62 2.0 235ms\n", " 13 -2.772170723015 + -14.21 -7.75 1.0 187ms\n", " 14 -2.772170723015 + -15.05 -8.16 1.0 214ms\n", "n Energy log10(ΔE) log10(Δρ) Diag Δtime\n", "--- --------------- --------- --------- ---- ------\n", " 1 -2.770625943749 -0.53 9.0 \n", " 2 -2.772045645794 -2.85 -1.30 1.0 171ms\n", " 3 -2.772082285558 -4.44 -2.70 1.0 188ms\n", " 4 -2.772083404636 -5.95 -3.52 2.0 191ms\n", " 5 -2.772083413612 -8.05 -3.74 2.0 191ms\n", " 6 -2.772083417735 -8.38 -5.31 1.0 217ms\n", " 7 -2.772083417809 -10.13 -5.40 3.0 219ms\n", " 8 -2.772083417811 -11.74 -5.94 1.0 178ms\n", " 9 -2.772083417811 -13.21 -6.34 1.0 204ms\n", " 10 -2.772083417811 -13.43 -6.86 1.0 181ms\n", " 11 -2.772083417811 + -13.88 -7.87 2.0 232ms\n", " 12 -2.772083417811 -14.12 -7.78 2.0 207ms\n", " 13 -2.772083417811 -14.07 -8.51 1.0 197ms\n" ] }, { "output_type": "execute_result", "data": { "text/plain": "1.773557979689211" }, "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.770789297183 -0.52 9.0 \n", " 2 -2.772057410798 -2.90 -1.32 1.0 166ms\n", " 3 -2.772083133636 -4.59 -2.52 1.0 171ms\n", " 4 -2.772083368767 -6.63 -3.32 1.0 210ms\n", " 5 -2.772083414757 -7.34 -3.78 2.0 195ms\n", " 6 -2.772083417668 -8.54 -4.67 1.0 174ms\n", " 7 -2.772083417798 -9.89 -4.98 2.0 259ms\n", " 8 -2.772083417810 -10.93 -5.68 1.0 182ms\n", " 9 -2.772083417811 -12.10 -6.09 2.0 216ms\n", " 10 -2.772083417811 -13.05 -7.21 1.0 206ms\n", " 11 -2.772083417811 + -14.21 -7.78 1.0 184ms\n", " 12 -2.772083417811 -14.03 -8.87 2.0 218ms\n", "\n", "Polarizability via ForwardDiff: 1.7725349677726303\n", "Polarizability via finite difference: 1.773557979689211\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.5" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.5", "language": "julia" } }, "nbformat": 4 }