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.08219381641361406Then 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.0006182370306135013
└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/terms/xc.jl:39
n Energy log10(ΔE) log10(Δρ) Magnet |Magn| Diag Δtime
--- --------------- --------- --------- ------ ------ ---- ------
1 -361.3854099976 0.07 1.334 3.440 6.9 3.78s
2 -363.2388880906 0.27 -0.21 0.014 3.623 3.4 3.65s
3 -363.3509262666 -0.95 -0.58 0.000 3.728 3.4 3.00s
4 -363.3889954954 -1.42 -1.18 0.000 3.717 2.6 2.46s
5 -363.3959769396 -2.16 -1.67 0.000 3.681 2.0 2.27s
6 -363.3973130826 -2.87 -2.04 0.000 3.657 1.5 3.21s
7 -363.3976107469 -3.53 -2.28 0.000 3.648 2.2 2.19s
8 -363.3976924516 -4.09 -2.64 0.000 3.647 1.6 2.02s
9 -363.3977070449 -4.84 -3.00 0.000 3.649 2.1 2.21s
10 -363.3977065659 + -6.32 -2.95 -0.000 3.649 2.0 2.04s
11 -363.3977091995 -5.58 -3.18 0.000 3.648 2.0 2.18s
12 -363.3977091772 + -7.65 -3.25 0.000 3.648 1.1 1.81s
13 -363.3977089780 + -6.70 -3.13 -0.000 3.648 2.0 2.11s
14 -363.3977086343 + -6.46 -3.03 -0.000 3.649 1.0 1.74s
15 -363.3977089399 -6.51 -2.94 -0.000 3.649 1.0 1.75s
16 -363.3977090564 -6.93 -2.94 -0.000 3.649 1.0 1.75s
17 -363.3977092316 -6.76 -2.95 0.000 3.649 1.0 1.75s
18 -363.3977093273 -7.02 -3.03 0.000 3.649 1.0 1.76s
19 -363.3977094677 -6.85 -3.11 0.000 3.649 1.0 1.76s
20 -363.3977094366 + -7.51 -3.06 0.000 3.649 1.0 1.76s
21 -363.3977095353 -7.01 -3.07 0.000 3.649 1.0 1.75s
22 -363.3977097969 -6.58 -3.12 0.000 3.649 1.0 1.75s
23 -363.3977098877 -7.04 -3.16 0.000 3.649 1.0 1.75s
24 -363.3977099013 -7.87 -3.30 0.000 3.649 1.0 1.76s
25 -363.3977100005 -7.00 -3.82 0.000 3.648 1.0 1.70s
26 -363.3977100137 -7.88 -4.07 0.000 3.648 1.2 1.74s
27 -363.3977100172 -8.45 -4.59 0.000 3.648 1.9 1.85s
28 -363.3977100178 -9.27 -4.95 0.000 3.648 2.8 2.13s
29 -363.3977100178 -10.22 -5.36 0.000 3.648 1.6 1.86s
30 -363.3977100178 -11.32 -5.34 0.000 3.648 1.6 3.15s
31 -363.3977100178 -11.08 -6.06 0.000 3.648 1.0 1.69sRun 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.11667630481949637With 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)