forked from JuliaMolSim/DFTK.jl
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathcompute_density.jl
136 lines (113 loc) · 6.16 KB
/
compute_density.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
using Test
using LinearAlgebra
using DFTK
import DFTK: total_local_potential, is_approx_integer
include("testcases.jl")
# TODO Once we have converged SCF densities in a file it would be better to instead / also
# test the energies of these densities and compare them directly to the reference
# energies obtained in the data files
if mpi_nprocs() == 1 # not easy to distribute
@testset "Using BZ symmetry yields identical density" begin
function get_bands(testcase, kgrid, kshift; symmetries, Ecut=5, tol=1e-8, n_rounds=1)
kwargs = ()
n_bands = div(testcase.n_electrons, 2, RoundUp)
if testcase.temperature !== nothing
kwargs = (temperature=testcase.temperature, smearing=DFTK.Smearing.FermiDirac())
n_bands = div(testcase.n_electrons, 2, RoundUp) + 4
end
occupation_threshold = 1e-7
model = model_DFT(testcase.lattice, testcase.atoms, testcase.positions,
:lda_xc_teter93; symmetries, kwargs...)
basis = PlaneWaveBasis(model; Ecut, kgrid, kshift, symmetries_respect_rgrid=false)
ham = Hamiltonian(basis; ρ=guess_density(basis))
res = diagonalize_all_kblocks(lobpcg_hyper, ham, n_bands; tol)
occ, εF = DFTK.compute_occupation(basis, res.λ; occupation_threshold)
ρnew = compute_density(basis, res.X, occ)
for it in 1:n_rounds
ham = Hamiltonian(basis; ρ=ρnew)
res = diagonalize_all_kblocks(lobpcg_hyper, ham, n_bands; tol=tol, ψguess=res.X)
occ, εF = DFTK.compute_occupation(basis, res.λ; occupation_threshold)
ρnew = compute_density(basis, res.X, occ)
end
ham, res.X, res.λ, ρnew, occ
end
function test_orthonormality(basis, ψ; tol=1e-8)
n_k = length(ψ)
n_states = size(ψ[1], 2)
n_fft = prod(basis.fft_size)
for (ik, kpt) in enumerate(basis.kpoints)
# Fourier-transform the wave functions to real space
ψk = ψ[ik]
ψk_real = cat((DFTK.G_to_r(basis, kpt, ψik) for ψik in eachcol(ψk))..., dims=4)
T = real(eltype(ψk_real))
ψk_real_mat = reshape(ψk_real, n_fft, n_states)
ψk_real_overlap = adjoint(ψk_real_mat) * ψk_real_mat
nondiag = ψk_real_overlap - I * (n_fft / basis.model.unit_cell_volume)
@test maximum(abs.(nondiag)) < tol
end
end
function test_full_vs_irreducible(testcase, kgrid_size; Ecut=5, tol=1e-8, n_ignore=0,
kshift=[0, 0, 0], eigenvectors=true)
res = get_bands(testcase, kgrid_size, kshift; symmetries=false, Ecut, tol)
ham_full, ψ_full, eigenvalues_full, ρ_full, occ_full = res
kfull = [kpt.coordinate for kpt in ham_full.basis.kpoints]
test_orthonormality(ham_full.basis, ψ_full; tol)
res = get_bands(testcase, kgrid_size, kshift; symmetries=true, Ecut, tol)
ham_ir, ψ_ir, eigenvalues_ir, ρ_ir, occ_ir = res
test_orthonormality(ham_ir.basis, ψ_ir; tol)
@test ham_full.basis.fft_size == ham_ir.basis.fft_size
# Test density is the same in both schemes, and symmetric wrt the basis symmetries
@test maximum(abs.(ρ_ir - ρ_full)) < 10tol
@test maximum(abs, DFTK.symmetrize_ρ(ham_ir.basis, ρ_ir) - ρ_ir) < tol
# Test local potential is the same in both schemes
@test maximum(abs, total_local_potential(ham_ir) - total_local_potential(ham_full)) < tol
# Test equivalent k-points have the same orbital energies
for (ik, kpt) in enumerate(ham_ir.basis.kpoints)
for symop in ham_ir.basis.symmetries
ikfull = findfirst(1:length(kfull)) do idx
all(is_approx_integer, kfull[idx] - symop.S * kpt.coordinate)
end
@test ikfull !== nothing
# Orbital energies should agree (by symmetry)
# The largest few are ignored, because the results between the
# full kpoints and the reduced kpoints are sometimes a little unstable
# (since we actually use different guesses)
range = 1:length(eigenvalues_full[ikfull]) - n_ignore
@test eigenvalues_full[ikfull][range] ≈ eigenvalues_ir[ik][range] atol=tol
end
end
if eigenvectors
# Test applying the symmetry transformation to the irreducible k-points
# yields an eigenfunction of the Hamiltonian
# Also check that the accumulated partial densities are equal
# to the returned density.
for (ik, kpt) in enumerate(ham_ir.basis.kpoints)
Hk_ir = ham_ir.blocks[ik]
for symop in ham_ir.basis.symmetries
Skpoint, ψSk = DFTK.apply_symop(symop, ham_ir.basis, Hk_ir.kpoint, ψ_ir[ik])
ikfull = findfirst(1:length(kfull)) do idx
all(is_approx_integer, round.(kfull[idx] - Skpoint.coordinate, digits=10))
end
@test !isnothing(ikfull)
Hk_full = ham_full.blocks[ikfull]
range = 1:length(eigenvalues_ir[ik]) - n_ignore
for iband in range
ψnSk = ψSk[:, iband]
residual = norm(Hk_full * ψnSk - eigenvalues_ir[ik][iband] * ψnSk)
@test residual < 10tol
end # iband
end # symop
end # k
end # eigenvectors
end
test_full_vs_irreducible(silicon, [3, 2, 3], Ecut=5, tol=1e-10)
test_full_vs_irreducible(silicon, [3, 3, 3], Ecut=5, tol=1e-6)
test_full_vs_irreducible(silicon, [3, 3, 3], Ecut=5, tol=1e-6, kshift=[1/2, 1/2, 1/2])
test_full_vs_irreducible(silicon, [3, 3, 3], Ecut=5, tol=1e-6, kshift=[1/2, 0, 1/2])
test_full_vs_irreducible(silicon, [3, 3, 3], Ecut=5, tol=1e-6, kshift=[0, 1/2, 0])
test_full_vs_irreducible(silicon, [2, 3, 4], Ecut=5, tol=1e-6)
test_full_vs_irreducible(magnesium, [2, 3, 4], Ecut=5, tol=1e-6, n_ignore=1)
test_full_vs_irreducible(aluminium, [1, 3, 5], Ecut=3, tol=1e-5, n_ignore=3,
eigenvectors=false)
end
end