22 use,
intrinsic :: iso_fortran_env, only: real64
47 integer :: n_simplices
48 integer,
allocatable :: simplices(:,:)
57 integer,
parameter :: submesh_tetras(6,4) = reshape([ &
63 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
81 function simplex_init(dim, naxis, nshifts, shift, kpoints, equiv)
result(this)
82 integer,
intent(in) :: dim
83 integer,
intent(in) :: naxis(1:dim)
84 integer,
intent(in) :: nshifts
85 real(real64),
intent(in) :: shift(:,:)
86 real(real64),
intent(in) :: kpoints(:,:)
87 integer,
intent(in) :: equiv(:)
88 type(simplex_t),
pointer :: this
90 real(real64) :: kmin(dim)
91 integer :: ik, npoints
94 integer,
allocatable :: kl123(:,:,:)
98 if (dim /= 3 .or. .not. all(naxis > 1))
then
99 message(1) =
"The linear tetrahedron method is currently only implemented for 3D systems"
102 if (nshifts /= 1)
then
103 message(1) =
"The linear tetrahedron method only works for automatic k-point grids with a single shift"
108 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
110 npoints = product(naxis)
111 kmin = minval(kpoints, 2)
114 ix(:) = nint((kpoints(:,ik) - kmin) * naxis + 1)
115 assert(kl123(ix(1), ix(2), ix(3)) == -1)
116 kl123(ix(1), ix(2), ix(3)) = equiv(ik)
119 this%n_simplices = 6 * npoints
121 safe_allocate(this%simplices(this%n_simplices, 4))
124 integer :: i, j, k, ip1, jp1, kp1, itetra, n
125 integer :: corners(8,3)
126 integer :: this_tetra(4), this_corner
130 ip1 = mod(i, naxis(1)) + 1
131 jp1 = mod(j, naxis(2)) + 1
132 kp1 = mod(k, naxis(3)) + 1
133 corners(:,:) = reshape([ &
141 ip1 , jp1 , kp1 ], shape(corners), order=[2, 1])
143 do itetra = 1,
size(submesh_tetras, 1)
144 n = (itetra - 1) + 6 * ((k - 1) + naxis(3) * ((j - 1) + naxis(2) * (i - 1))) + 1
145 this_tetra(:) = submesh_tetras(itetra, :)
147 this_corner = kl123(corners(this_tetra(ik), 1), corners(this_tetra(ik), 2), corners(this_tetra(ik), 3))
148 this%simplices(n,ik) = this_corner
156 safe_deallocate_a(kl123)
167 safe_deallocate_a(this%simplices)
181 real(real64),
intent(in) :: etetra(4)
182 real(real64),
intent(in) :: ef
183 real(real64),
intent(out) :: weights(4)
184 real(real64),
intent(out) :: dos
186 real(real64) :: e(4), e1, e2, e3, e4
190 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
191 real(real64),
parameter :: vt_4vg = vt_vg / 4.0_real64
206 elseif (e4 < ef)
then
209 elseif (e1 < ef .and. ef <= e2)
then
211 real(real64) :: e21, e31, e41, c
215 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
218 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
223 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
225 elseif (e2 < ef .and. ef <= e3)
then
227 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
233 c1 = vt_4vg * (ef - e1) ** 2 / (e41 * e31)
234 c2 = vt_4vg * (ef - e1) * (ef - e2) * (e3 - ef) / (e41 * e32 * e31)
235 c3 = vt_4vg * (ef - e2) ** 2 * (e4 - ef) / (e42 * e32 * e41)
238 c1 + (c1 + c2) * (e3 - ef) / e31 + (c1 + c2 + c3) * (e4 - ef) / e41, &
239 c1 + c2 + c3 + (c2 + c3) * (e3 - ef) / e32 + c3 * (e4 - ef) / e42, &
240 (c1 + c2) * (ef - e1) / e31 + (c2 + c3) * (ef - e2) / e32, &
241 (c1 + c2 + c3) * (ef - e1) / e41 + c3 * (ef - e2) / e42]
243 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
244 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
246 elseif (e3 < ef .and. ef <= e4)
then
248 real(real64) :: e41, e42, e43, c
252 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
254 w(:) = vt_4vg - c * [ &
258 4.0_real64 - (e4 - ef) * (m_one / e41 + m_one / e42 + m_one / e43)]
260 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
266 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
278 real(real64),
intent(in) :: etetra(4)
279 real(real64),
intent(in) :: ef
280 real(real64),
intent(out) :: dos
282 real(real64) :: e(4), e1, e2, e3, e4
284 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
295 if (e1 >= ef .or. e4 < ef)
then
297 elseif (e1 < ef .and. ef <= e2)
then
299 real(real64) :: e21, e31, e41
303 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
305 elseif (e2 < ef .and. ef <= e3)
then
307 real(real64) :: e21, e31, e32, e41, e42
313 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
314 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
316 elseif (e3 < ef .and. ef <= e4)
then
318 real(real64) :: e41, e42, e43
322 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
This is the common interface to a sorting routine. It performs the shell algorithm,...
type(debug_t), save, public debug
subroutine, public debug_pop_sub(sub_name)
Pop a routine from the debug trace.
subroutine, public debug_push_sub(sub_name)
Push a routine to the debug trace.
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
integer(int64), public global_sizeof
logical pure function, public not_in_openmp()
character(len=100), public global_alloc_errmsg
integer, public global_alloc_err
real(real64), parameter, public m_one
subroutine, public alloc_error(size, file, line)
subroutine, public dealloc_error(size, file, line)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
type(profile_vars_t), target, save, public prof_vars
integer, parameter, public profiling_memory
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine, public profiling_memory_deallocate(var, file, line, size)
subroutine, public profiling_memory_allocate(var, file, line, size_)
subroutine, public simplex_dos_3d(etetra, eF, dos)
Get only the DOS contribution of a single tetrahedron.
subroutine, public simplex_end(this)
Destructor for linear simplex methods.
type(simplex_t) function, pointer, public simplex_init(dim, naxis, nshifts, shift, kpoints, equiv)
Constructor for linear simplex methods.
subroutine, public simplex_weights_3d(etetra, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
This module is intended to contain "only mathematical" functions and procedures.