{ "cells": [ { "cell_type": "markdown", "source": [ "# Polarizability using automatic differentiation\n", "\n", "Simple example for computing properties using (forward-mode)\n", "automatic differentation.\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=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.770907117738 -0.52 9.0\n", " 2 -2.772149038433 -2.91 -1.32 1.0\n", " 3 -2.772170030240 -4.68 -2.39 2.0\n", " 4 -2.772170722605 -6.16 -4.07 2.0\n", "┌ Warning: Negative ρ detected\n", "│ min_ρ = -1.237664896679902e-18\n", "└ @ DFTK /home/runner/work/DFTK.jl/DFTK.jl/src/densities.jl:7\n", " 5 -2.772170723005 -9.40 -5.00 2.0\n", "n Energy log10(ΔE) log10(Δρ) Diag\n", "--- --------------- --------- --------- ----\n", " 1 -2.770681324509 -0.53 8.0\n", " 2 -2.772036828676 -2.87 -1.30 1.0\n", " 3 -2.772083327254 -4.33 -2.78 2.0\n", " 4 -2.772083416680 -7.05 -3.80 2.0\n", " 5 -2.772083417796 -8.95 -4.60 2.0\n" ] }, { "output_type": "execute_result", "data": { "text/plain": "1.77364460650267" }, "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.770741463218 -0.53 8.0\n", " 2 -2.772049843356 -2.88 -1.30 1.0\n", " 3 -2.772082720625 -4.48 -2.69 1.0\n", " 4 -2.772083413349 -6.16 -3.77 2.0\n", " 5 -2.772083417381 -8.39 -4.22 2.0\n", "\n", "Polarizability via ForwardDiff: 1.7724740060017201\n", "Polarizability via finite difference: 1.77364460650267\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.7.2" }, "kernelspec": { "name": "julia-1.7", "display_name": "Julia 1.7.2", "language": "julia" } }, "nbformat": 4 }