22 use,
intrinsic :: iso_fortran_env, only: real64
49 integer :: n_simplices
50 integer,
allocatable :: simplices(:,:)
75 function simplex_init(dim, naxis, nshifts, shift, kpoints, equiv)
result(this)
76 integer,
intent(in) :: dim
77 integer,
intent(in) :: naxis(1:dim)
78 integer,
intent(in) :: nshifts
79 real(real64),
intent(in) :: shift(:,:)
80 real(real64),
intent(in) :: kpoints(:,:)
81 integer,
intent(in) :: equiv(:)
82 type(simplex_t),
pointer :: this
84 real(real64) :: kmin(dim)
85 integer :: ik, npoints
88 integer,
allocatable :: kl123(:,:,:)
90 integer :: rdim, raxis(3)
94 if (nshifts /= 1)
then
95 message(1) =
"The linear tetrahedron method only works for automatic k-point grids with a single shift"
100 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
102 npoints = product(naxis)
103 kmin = minval(kpoints, 2)
106 ix(:) = nint((kpoints(:,ik) - kmin) * naxis + 1)
107 assert(kl123(ix(1), ix(2), ix(3)) == -1)
108 kl123(ix(1), ix(2), ix(3)) = equiv(ik)
111 rdim = sum(merge(1, 0, naxis > 1))
112 raxis(1:rdim) = pack(naxis, naxis > 1)
114 if (any(raxis(1:rdim) /= naxis(1:rdim)))
then
115 message(1) =
"The periodic dimensions must be consecutive"
126 integer,
parameter :: submesh_triangles(2,3) = reshape([ &
128 1, 4, 3], shape(submesh_triangles), order=[2, 1])
130 integer :: i, j, ip1, jp1, it, n
131 integer :: corners(4,2)
132 integer :: this_triangle(3), this_corner
134 this%n_simplices = 2 * npoints
136 safe_allocate(this%simplices(this%n_simplices, this%dim))
140 ip1 = mod(i, raxis(1)) + 1
141 jp1 = mod(j, raxis(2)) + 1
142 corners(:,:) = reshape([ &
146 i , jp1 ], shape(corners), order=[2, 1])
148 do it = 1,
size(submesh_triangles, 1)
149 n = (it - 1) + 2 * ((j - 1) + raxis(2) * (i - 1)) + 1
150 this_triangle(:) = submesh_triangles(it, :)
152 this_corner = kl123(corners(this_triangle(ik), 1), corners(this_triangle(ik), 2), 1)
153 this%simplices(n,ik) = this_corner
161 integer,
parameter :: submesh_tetras(6,4) = reshape([ &
167 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
169 integer :: i, j, k, ip1, jp1, kp1, it, n
170 integer :: corners(8,3)
171 integer :: this_tetra(4), this_corner
173 this%n_simplices = 6 * npoints
175 safe_allocate(this%simplices(this%n_simplices, this%dim))
180 ip1 = mod(i, raxis(1)) + 1
181 jp1 = mod(j, raxis(2)) + 1
182 kp1 = mod(k, raxis(3)) + 1
183 corners(:,:) = reshape([ &
191 ip1 , jp1 , kp1 ], shape(corners), order=[2, 1])
193 do it = 1,
size(submesh_tetras, 1)
194 n = (it - 1) + 6 * ((k - 1) + raxis(3) * ((j - 1) + raxis(2) * (i - 1))) + 1
195 this_tetra(:) = submesh_tetras(it, :)
197 this_corner = kl123(corners(this_tetra(ik), 1), corners(this_tetra(ik), 2), corners(this_tetra(ik), 3))
198 this%simplices(n,ik) = this_corner
207 safe_deallocate_a(kl123)
218 safe_deallocate_a(this%simplices)
232 real(real64),
intent(in) :: etriangle(3)
233 real(real64),
intent(in) :: ef
234 real(real64),
intent(out) :: weights(3)
235 real(real64),
intent(out) :: dos
237 real(real64) :: e(3), e1, e2, e3
241 real(real64),
parameter :: vt_vg = 1.0_real64 / 2.0_real64
242 real(real64),
parameter :: vt_3vg = vt_vg / 3.0_real64
256 elseif (e3 < ef)
then
259 elseif (e1 < ef .and. ef <= e2)
then
261 real(real64) :: e21, e31, c
264 c = vt_3vg * (ef - e1) ** 2 / (e21 * e31)
267 3.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31), &
271 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
273 elseif (e2 < ef .and. ef <= e3)
then
275 real(real64) :: e23, e31, c1, c2
279 c2 = vt_3vg * (ef - e3) ** 2 / (e23 * e31)
282 c1 - c2 * (ef - e3) / e31, &
283 c1 + c2 * (ef - e3) / e23, &
284 c1 + c2 * (3.0_real64 - (ef - e3) * (m_one / e23 - m_one / e31))]
286 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
292 weights(idx) = w + dos * (sum(e) - 3.0_real64 * e) / 24.0_real64
301 real(real64),
intent(in) :: etriangle(3)
302 real(real64),
intent(in) :: ef
303 real(real64),
intent(out) :: dos
305 real(real64) :: e(3), e1, e2, e3
307 real(real64),
parameter :: vt_vg = 1.0_real64 / 2.0_real64
317 if (ef <= e1 .or. e3 < ef)
then
319 elseif (e1 < ef .and. ef <= e2)
then
321 real(real64) :: e21, e31
324 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
326 elseif (e2 < ef .and. ef <= e3)
then
328 real(real64) :: e23, e31
331 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
348 real(real64),
intent(in) :: etetra(4)
349 real(real64),
intent(in) :: ef
350 real(real64),
intent(out) :: weights(4)
351 real(real64),
intent(out) :: dos
353 real(real64) :: e(4), e1, e2, e3, e4
357 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
358 real(real64),
parameter :: vt_4vg = vt_vg / 4.0_real64
373 elseif (e4 < ef)
then
376 elseif (e1 < ef .and. ef <= e2)
then
378 real(real64) :: e21, e31, e41, c
382 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
385 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
390 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
392 elseif (e2 < ef .and. ef <= e3)
then
394 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
400 c1 = vt_4vg * (ef - e1) ** 2 / (e41 * e31)
401 c2 = vt_4vg * (ef - e1) * (ef - e2) * (e3 - ef) / (e41 * e32 * e31)
402 c3 = vt_4vg * (ef - e2) ** 2 * (e4 - ef) / (e42 * e32 * e41)
405 c1 + (c1 + c2) * (e3 - ef) / e31 + (c1 + c2 + c3) * (e4 - ef) / e41, &
406 c1 + c2 + c3 + (c2 + c3) * (e3 - ef) / e32 + c3 * (e4 - ef) / e42, &
407 (c1 + c2) * (ef - e1) / e31 + (c2 + c3) * (ef - e2) / e32, &
408 (c1 + c2 + c3) * (ef - e1) / e41 + c3 * (ef - e2) / e42]
410 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
411 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
413 elseif (e3 < ef .and. ef <= e4)
then
415 real(real64) :: e41, e42, e43, c
419 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
421 w(:) = vt_4vg - c * [ &
425 4.0_real64 - (e4 - ef) * (m_one / e41 + m_one / e42 + m_one / e43)]
427 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
433 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
445 real(real64),
intent(in) :: etetra(4)
446 real(real64),
intent(in) :: ef
447 real(real64),
intent(out) :: dos
449 real(real64) :: e(4), e1, e2, e3, e4
451 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
462 if (e1 >= ef .or. e4 < ef)
then
464 elseif (e1 < ef .and. ef <= e2)
then
466 real(real64) :: e21, e31, e41
470 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
472 elseif (e2 < ef .and. ef <= e3)
then
474 real(real64) :: e21, e31, e32, e41, e42
480 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
481 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
483 elseif (e3 < ef .and. ef <= e4)
then
485 real(real64) :: e41, e42, e43
489 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 messages_not_implemented(feature, namespace)
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_weights_2d(etriangle, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
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_dos_2d(etriangle, eF, dos)
Get only the DOS contribution of a single triangle.
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.