Hubbard correction (DFT+U)
In this example, we'll plot the DOS and projected DOS of Nickel Oxide with and without the Hubbard term correction.
using DFTK
using PseudoPotentialData
using Unitful
using UnitfulAtomic
using PlotsDefine the geometry and pseudopotential
a = 7.9 # Nickel Oxide lattice constant in Bohr
lattice = a * [[ 1.0 0.5 0.5];
[ 0.5 1.0 0.5];
[ 0.5 0.5 1.0]]
pseudopotentials = PseudoFamily("dojo.nc.sr.pbe.v0_4_1.standard.upf")
Ni = ElementPsp(:Ni, pseudopotentials)
O = ElementPsp(:O, pseudopotentials)
atoms = [Ni, O, Ni, O]
positions = [zeros(3), ones(3) / 4, ones(3) / 2, ones(3) * 3 / 4]
magnetic_moments = [2, 0, -1, 0]4-element Vector{Int64}:
2
0
-1
0First, we run an SCF and band computation without the Hubbard term
model = model_DFT(lattice, atoms, positions; temperature=5e-3,
functionals=PBE(), magnetic_moments)
basis = PlaneWaveBasis(model; Ecut=20, kgrid=[2, 2, 2])
scfres = self_consistent_field(basis; tol=1e-6, ρ=guess_density(basis, magnetic_moments))
bands = compute_bands(scfres, MonkhorstPack(4, 4, 4))
lowest_unocc_band = findfirst(ε -> ε-bands.εF > 0, bands.eigenvalues[1])
band_gap = bands.eigenvalues[1][lowest_unocc_band] - bands.eigenvalues[1][lowest_unocc_band-1]0.08219297178319335Then we plot the DOS and the PDOS for the relevant 3D (pseudo)atomic projector
εF = bands.εF
width = 5.0u"eV"
εrange = (εF - austrip(width), εF + austrip(width))
p = plot_dos(bands; εrange, colors=[:red, :red])
plot_pdos(bands; p, iatom=1, label="3D", colors=[:yellow, :orange], εrange)To perform and Hubbard computation, we have to define the Hubbard manifold and associated constant.
In DFTK there are a few ways to construct the OrbitalManifold. Here, we will apply the Hubbard correction on the 3D orbital of all nickel atoms. To select all nickel atoms, we can:
- Pass the
Nielement directly. - Pass the
:Nisymbol. - Pass the list of atom indices, here
[1, 3].
To select the orbitals, it is recommended to use their label, such as "3D" for PseudoDojo pseudopotentials.
Note that "manifold" is the standard term used in the literature for the set of atomic orbitals used to compute the Hubbard correction, but it is not meant in the mathematical sense.
U = 10u"eV"
# Alternative:
# manifold = OrbitalManifold(:Ni, "3D")
# Alternative:
# manifold = OrbitalManifold([1, 3], "3D")
manifold = OrbitalManifold(Ni, "3D")OrbitalManifold(Ni, "3D")Run SCF with a DFT+U setup, notice the extra_terms keyword argument, setting up the Hubbard +U term.
model = model_DFT(lattice, atoms, positions; extra_terms=[Hubbard(manifold, U)],
functionals=PBE(), temperature=5e-3, magnetic_moments)
basis = PlaneWaveBasis(model; Ecut=20, kgrid=[2, 2, 2])
scfres = self_consistent_field(basis; tol=1e-6, ρ=guess_density(basis, magnetic_moments));┌ Warning: Negative ρcore detected: -0.0006182370306135084
└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/terms/xc.jl:39
n Energy log10(ΔE) log10(Δρ) Magnet |Magn| Diag Δtime
--- --------------- --------- --------- ------ ------ ---- ------
1 -361.3863537750 0.07 1.334 3.441 6.8 3.68s
2 -363.2379453148 0.27 -0.21 0.014 3.624 3.2 6.02s
3 -363.3509555629 -0.95 -0.58 0.000 3.727 3.1 2.90s
4 -363.3890101251 -1.42 -1.18 0.000 3.717 2.6 2.45s
5 -363.3959857415 -2.16 -1.67 0.000 3.681 2.0 2.27s
6 -363.3973057607 -2.88 -2.04 0.000 3.656 1.5 1.91s
7 -363.3976115193 -3.51 -2.28 0.000 3.647 2.4 2.21s
8 -363.3976872573 -4.12 -2.57 0.000 3.647 1.5 1.96s
9 -363.3977062299 -4.72 -2.97 0.000 3.649 2.0 2.15s
10 -363.3977061309 + -7.00 -2.91 -0.000 3.649 2.6 2.09s
11 -363.3977087968 -5.57 -3.11 -0.000 3.649 1.2 1.82s
12 -363.3977090930 -6.53 -3.03 -0.000 3.648 1.0 1.79s
13 -363.3977096388 -6.26 -2.97 0.000 3.648 1.4 1.83s
14 -363.3977086280 + -6.00 -2.86 0.000 3.648 1.9 1.98s
15 -363.3977093680 -6.13 -2.88 0.000 3.648 1.0 1.80s
16 -363.3977095514 -6.74 -2.87 -0.000 3.648 1.0 1.78s
17 -363.3977094982 + -7.27 -2.85 0.000 3.648 1.0 1.78s
18 -363.3977094842 + -7.85 -2.83 0.000 3.648 1.0 1.78s
19 -363.3977094313 + -7.28 -2.82 0.000 3.648 1.0 1.78s
20 -363.3977093200 + -6.95 -2.81 0.000 3.648 1.0 1.79s
21 -363.3977092221 + -7.01 -2.81 0.000 3.648 1.0 1.78s
22 -363.3977091710 + -7.29 -2.80 0.000 3.648 1.0 1.78s
23 -363.3977095548 -6.42 -2.63 -0.000 3.648 1.0 1.77s
24 -363.3977095310 + -7.62 -2.58 -0.000 3.648 1.0 1.78s
25 -363.3977095933 -7.21 -2.73 -0.000 3.648 1.0 1.78s
26 -363.3977099612 -6.43 -3.32 -0.000 3.648 1.0 3.00s
27 -363.3977099375 + -7.62 -3.58 -0.000 3.649 1.0 1.64s
28 -363.3977099938 -7.25 -3.27 -0.000 3.648 1.0 1.66s
29 -363.3977100141 -7.69 -4.16 -0.000 3.648 1.0 1.65s
30 -363.3977100064 + -8.11 -4.12 -0.000 3.648 1.9 1.92s
31 -363.3977100136 -8.14 -4.29 -0.000 3.648 1.0 1.68s
32 -363.3977100164 -8.55 -4.29 -0.000 3.648 1.4 1.72s
33 -363.3977100167 -9.55 -4.55 -0.000 3.648 1.0 1.67s
34 -363.3977100171 -9.39 -4.78 -0.000 3.648 1.2 1.70s
35 -363.3977100177 -9.25 -5.10 0.000 3.648 1.1 1.69s
36 -363.3977100178 -9.80 -5.36 0.000 3.648 2.0 1.96s
37 -363.3977100178 -10.63 -5.67 0.000 3.648 2.0 1.90s
38 -363.3977100179 -11.16 -6.18 0.000 3.648 1.0 1.67sRun band computation
bands_hub = compute_bands(scfres, MonkhorstPack(4, 4, 4))
lowest_unocc_band = findfirst(ε -> ε-bands_hub.εF > 0, bands_hub.eigenvalues[1])
band_gap = bands_hub.eigenvalues[1][lowest_unocc_band] - bands_hub.eigenvalues[1][lowest_unocc_band-1]0.1166761992968941With the electron localization introduced by the Hubbard term, the band gap has now opened, reflecting the experimental insulating behaviour of Nickel Oxide.
εF = bands_hub.εF
εrange = (εF - austrip(width), εF + austrip(width))
p = plot_dos(bands_hub; p, colors=[:blue, :blue], εrange)
plot_pdos(bands_hub; p, iatom=1, label="3D", colors=[:green, :purple], εrange)