22 use,
intrinsic :: iso_fortran_env, only: real64
49 integer :: n_simplices
50 integer,
allocatable :: simplices(:,:)
76 function simplex_init(dim, naxis, nshifts, shift, kpoints, equiv, opt)
result(this)
77 integer,
intent(in) :: dim
78 integer,
intent(in) :: naxis(1:dim)
79 integer,
intent(in) :: nshifts
80 real(real64),
intent(in) :: shift(:,:)
81 real(real64),
intent(in) :: kpoints(:,:)
82 integer,
intent(in) :: equiv(:)
83 logical,
intent(in) :: opt
84 type(simplex_t),
pointer :: this
86 real(real64) :: kmin(dim)
87 integer :: ik, npoints
90 integer,
allocatable :: kl123(:,:,:)
92 integer :: rdim, raxis(3)
96 if (nshifts /= 1)
then
97 message(1) =
"The linear tetrahedron method only works for automatic k-point grids with a single shift"
102 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
104 npoints = product(naxis)
105 kmin = minval(kpoints, 2)
108 ix(:) = nint((kpoints(:,ik) - kmin) * naxis + 1)
109 assert(kl123(ix(1), ix(2), ix(3)) == -1)
110 kl123(ix(1), ix(2), ix(3)) = equiv(ik)
113 rdim = sum(merge(1, 0, naxis > 1))
114 raxis(1:rdim) = pack(naxis, naxis > 1)
116 if (any(raxis(1:rdim) /= naxis(1:rdim)))
then
117 message(1) =
"The periodic dimensions must be consecutive"
128 integer,
parameter :: submesh_triangles(2,3) = reshape([ &
130 1, 4, 3], shape(submesh_triangles), order=[2, 1])
132 integer,
parameter :: b(10,3) = reshape([ &
143 ], shape(b), order=[2, 1])
145 integer :: i, j, ip1, jp1, it, n
146 integer :: corners(4,2), v(3,2), c(2)
147 integer :: this_triangle(3), this_corner
149 this%n_simplices = 2 * npoints
150 this%dim = merge(10, 3, opt)
151 safe_allocate(this%simplices(this%n_simplices, this%dim))
155 ip1 = modulo(i, raxis(1)) + 1
156 jp1 = modulo(j, raxis(2)) + 1
157 corners(:,:) = reshape([ &
161 i , jp1 ], shape(corners), order=[2, 1])
163 do it = 1,
size(submesh_triangles, 1)
164 n = (it - 1) + 2 * ((j - 1) + raxis(2) * (i - 1)) + 1
165 this_triangle(:) = submesh_triangles(it, :)
166 v(1,:) = corners(this_triangle(1), :)
167 v(2,:) = corners(this_triangle(2), :)
168 v(3,:) = corners(this_triangle(3), :)
170 c(:) = b(ik,1) * v(1,:) + b(ik,2) * v(2,:) + b(ik,3) * v(3,:)
171 c(:) = modulo(c(:) - 1, raxis(1:rdim)) + 1
172 this_corner = kl123(c(1), c(2), 1)
173 this%simplices(n,ik) = this_corner
181 integer,
parameter :: submesh_tetras(6,4) = reshape([ &
187 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
189 integer,
parameter :: b(20,4) = reshape([ &
210 ], shape(b), order=[2, 1])
212 integer :: i, j, k, ip1, jp1, kp1, it, n
213 integer :: corners(8,3), v(4,3), c(3)
214 integer :: this_tetra(4), this_corner
216 this%n_simplices = 6 * npoints
217 this%dim = merge(20, 4, opt)
218 safe_allocate(this%simplices(this%n_simplices, this%dim))
223 ip1 = modulo(i, raxis(1)) + 1
224 jp1 = modulo(j, raxis(2)) + 1
225 kp1 = modulo(k, raxis(3)) + 1
226 corners(:,:) = reshape([ &
234 ip1 , jp1 , kp1 ], shape(corners), order=[2, 1])
236 do it = 1,
size(submesh_tetras, 1)
237 n = (it - 1) + 6 * ((k - 1) + raxis(3) * ((j - 1) + raxis(2) * (i - 1))) + 1
238 this_tetra(:) = submesh_tetras(it, :)
239 v(1,:) = corners(this_tetra(1), :)
240 v(2,:) = corners(this_tetra(2), :)
241 v(3,:) = corners(this_tetra(3), :)
242 v(4,:) = corners(this_tetra(4), :)
244 c(:) = b(ik,1) * v(1,:) + b(ik,2) * v(2,:) + b(ik,3) * v(3,:) + b(ik,4) * v(4,:)
245 c(:) = modulo(c(:) - 1, raxis(1:rdim)) + 1
246 this_corner = kl123(c(1), c(2), c(3))
247 this%simplices(n,ik) = this_corner
256 safe_deallocate_a(kl123)
267 safe_deallocate_a(this%simplices)
281 real(real64),
intent(in) :: etriangle(:)
282 real(real64),
intent(in) :: ef
283 real(real64),
intent(out) :: weights(3)
284 real(real64),
intent(out) :: dos
286 real(real64) :: e(3), e1, e2, e3
290 real(real64),
parameter :: vt_vg = 1.0_real64 / 2.0_real64
291 real(real64),
parameter :: vt_3vg = vt_vg / 3.0_real64
293 real(real64),
parameter :: p(3,10) = 1.0_real64 / 360.0_real64 * reshape([ &
294 402 , 0 , 6 , -13 , 5 , -17 , -13 , -11 , 7 , -6 , &
295 6 , 396 , 6 , -9 , -15 , 3 , 3 , -15 , -9 , -6 , &
296 6 , 0 , 402 , 7 , -11 , -13 , -17 , 5 , -13 , -6 &
297 ], shape(p), order=[2, 1])
301 select case (
size(etriangle))
308 do i = 1,
size(etriangle)
309 e(:) = e(:) + p(:,i) * etriangle(i)
325 elseif (e3 < ef)
then
328 elseif (e1 < ef .and. ef <= e2)
then
330 real(real64) :: e21, e31, c
333 c = vt_3vg * (ef - e1) ** 2 / (e21 * e31)
336 3.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31), &
340 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
342 elseif (e2 < ef .and. ef <= e3)
then
344 real(real64) :: e23, e31, c1, c2
348 c2 = vt_3vg * (ef - e3) ** 2 / (e23 * e31)
351 c1 - c2 * (ef - e3) / e31, &
352 c1 + c2 * (ef - e3) / e23, &
353 c1 + c2 * (3.0_real64 - (ef - e3) * (m_one / e23 - m_one / e31))]
355 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
361 select case (
size(etriangle))
363 weights(idx) = w + dos * (sum(e) - 3.0_real64 * e) / 24.0_real64
377 real(real64),
intent(in) :: etriangle(:)
378 real(real64),
intent(in) :: ef
379 real(real64),
intent(out) :: dos
381 real(real64) :: e(3), e1, e2, e3
383 real(real64),
parameter :: vt_vg = 1.0_real64 / 2.0_real64
385 real(real64),
parameter :: p(3,10) = 1.0_real64 / 360.0_real64 * reshape([ &
386 402 , 0 , 6 , -13 , 5 , -17 , -13 , -11 , 7 , -6 , &
387 6 , 396 , 6 , -9 , -15 , 3 , 3 , -15 , -9 , -6 , &
388 6 , 0 , 402 , 7 , -11 , -13 , -17 , 5 , -13 , -6 &
389 ], shape(p), order=[2, 1])
393 select case (
size(etriangle))
400 do i = 1,
size(etriangle)
401 e(:) = e(:) + p(:,i) * etriangle(i)
413 if (ef <= e1 .or. e3 < ef)
then
415 elseif (e1 < ef .and. ef <= e2)
then
417 real(real64) :: e21, e31
420 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
422 elseif (e2 < ef .and. ef <= e3)
then
424 real(real64) :: e23, e31
427 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
447 real(real64),
intent(in) :: etetra(:)
448 real(real64),
intent(in) :: ef
449 real(real64),
intent(out) :: weights(4)
450 real(real64),
intent(out) :: dos
452 real(real64) :: e(4), e1, e2, e3, e4
456 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
457 real(real64),
parameter :: vt_4vg = vt_vg / 4.0_real64
459 real(real64),
parameter :: p(4,20) = 1.0_real64 / 1260.0_real64 * reshape([ &
460 1440, 0 , 30 , 0 , -38 , 7 , 17 , -28 , -56 , 9 , -46 , 9 , -38 , -28 , 17 , 7 , -18 , -18 , 12 , -18 , &
461 0 , 1440, 0 , 30 , -28 , -38 , 7 , 17 , 9 , -56 , 9 , -46 , 7 , -38 , -28 , 17 , -18 , -18 , -18 , 12 , &
462 30 , 0 , 1440, 0 , 17 , -28 , -38 , 7 , -46 , 9 , -56 , 9 , 17 , 7 , -38 , -28 , 12 , -18 , -18 , -18 , &
463 0 , 30 , 0 , 1440, 7 , 17 , -28 , -38 , 9 , -46 , 9 , -56 , -28 , 17 , 7 , -38 , -18 , 12 , -18 , -18 &
464 ], shape(p), order=[2, 1])
468 select case (
size(etetra))
475 do i = 1,
size(etetra)
476 e(:) = e(:) + p(:,i) * etetra(i)
493 elseif (e4 < ef)
then
496 elseif (e1 < ef .and. ef <= e2)
then
498 real(real64) :: e21, e31, e41, c
502 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
505 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
510 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
512 elseif (e2 < ef .and. ef <= e3)
then
514 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
520 c1 = vt_4vg * (ef - e1) ** 2 / (e41 * e31)
521 c2 = vt_4vg * (ef - e1) * (ef - e2) * (e3 - ef) / (e41 * e32 * e31)
522 c3 = vt_4vg * (ef - e2) ** 2 * (e4 - ef) / (e42 * e32 * e41)
525 c1 + (c1 + c2) * (e3 - ef) / e31 + (c1 + c2 + c3) * (e4 - ef) / e41, &
526 c1 + c2 + c3 + (c2 + c3) * (e3 - ef) / e32 + c3 * (e4 - ef) / e42, &
527 (c1 + c2) * (ef - e1) / e31 + (c2 + c3) * (ef - e2) / e32, &
528 (c1 + c2 + c3) * (ef - e1) / e41 + c3 * (ef - e2) / e42]
530 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
531 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
533 elseif (e3 < ef .and. ef <= e4)
then
535 real(real64) :: e41, e42, e43, c
539 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
541 w(:) = vt_4vg - c * [ &
545 4.0_real64 - (e4 - ef) * (m_one / e41 + m_one / e42 + m_one / e43)]
547 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
553 select case (
size(etetra))
555 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
575 real(real64),
intent(in) :: etetra(:)
576 real(real64),
intent(in) :: ef
577 real(real64),
intent(out) :: dos
579 real(real64) :: e(4), e1, e2, e3, e4
581 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
583 real(real64),
parameter :: p(4,20) = 1.0_real64 / 1260.0_real64 * reshape([ &
584 1440, 0 , 30 , 0 , -38 , 7 , 17 , -28 , -56 , 9 , -46 , 9 , -38 , -28 , 17 , 7 , -18 , -18 , 12 , -18 , &
585 0 , 1440, 0 , 30 , -28 , -38 , 7 , 17 , 9 , -56 , 9 , -46 , 7 , -38 , -28 , 17 , -18 , -18 , -18 , 12 , &
586 30 , 0 , 1440, 0 , 17 , -28 , -38 , 7 , -46 , 9 , -56 , 9 , 17 , 7 , -38 , -28 , 12 , -18 , -18 , -18 , &
587 0 , 30 , 0 , 1440, 7 , 17 , -28 , -38 , 9 , -46 , 9 , -56 , -28 , 17 , 7 , -38 , -18 , 12 , -18 , -18 &
588 ], shape(p), order=[2, 1])
592 select case (
size(etetra))
599 do i = 1,
size(etetra)
600 e(:) = e(:) + p(:,i) * etetra(i)
613 if (e1 >= ef .or. e4 < ef)
then
615 elseif (e1 < ef .and. ef <= e2)
then
617 real(real64) :: e21, e31, e41
621 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
623 elseif (e2 < ef .and. ef <= e3)
then
625 real(real64) :: e21, e31, e32, e41, e42
631 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
632 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
634 elseif (e3 < ef .and. ef <= e4)
then
636 real(real64) :: e41, e42, e43
640 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.
type(simplex_t) function, pointer, public simplex_init(dim, naxis, nshifts, shift, kpoints, equiv, opt)
Constructor for linear simplex methods.
subroutine, public simplex_end(this)
Destructor 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.