22 use,
intrinsic :: iso_fortran_env, only: real64
54 integer :: n_simplices
55 integer,
allocatable :: simplices(:,:)
81 function simplex_init(dim, naxis, nshifts, shift, kpoints, equiv, opt)
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 logical,
intent(in) :: opt
89 type(simplex_t),
pointer :: this
91 real(real64) :: kmin(dim)
92 integer :: ik, npoints
95 integer,
allocatable :: kl123(:,:,:)
97 integer :: rdim, raxis(3)
101 if (nshifts /= 1)
then
102 message(1) =
"The linear tetrahedron method only works for automatic k-point grids with a single shift"
107 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
109 npoints = product(naxis)
110 kmin = minval(kpoints, 2)
113 ix(:) = nint((kpoints(:,ik) - kmin) * naxis + 1)
114 assert(kl123(ix(1), ix(2), ix(3)) == -1)
115 kl123(ix(1), ix(2), ix(3)) = equiv(ik)
118 rdim = sum(merge(1, 0, naxis > 1))
119 raxis(1:rdim) = pack(naxis, naxis > 1)
121 if (any(raxis(1:rdim) /= naxis(1:rdim)))
then
122 message(1) =
"The periodic dimensions must be consecutive"
131 integer,
parameter :: submesh_segments(1,2) = reshape([ &
132 1, 2 ], shape(submesh_segments), order=[2, 1])
134 integer,
parameter :: b(4,2) = reshape([ &
139 ], shape(b), order=[2, 1])
141 integer :: i, ip1, it, n
142 integer :: corners(2,1), v(2,1), c(1)
143 integer :: this_segment(2), this_corner
145 this%n_simplices = npoints
147 this%sdim = merge(4, 2, opt)
148 safe_allocate(this%simplices(this%n_simplices, this%sdim))
151 ip1 = modulo(i, raxis(1)) + 1
152 corners(:,:) = reshape([ i , ip1 ], shape(corners), order=[2, 1])
154 do it = 1,
size(submesh_segments, 1)
155 n = (it - 1) + 1 * (i - 1) + 1
156 this_segment(:) = submesh_segments(it, :)
157 v(1,:) = corners(this_segment(1), :)
158 v(2,:) = corners(this_segment(2), :)
160 c(:) = b(ik,1) * v(1,:) + b(ik,2) * v(2,:)
161 c(:) = modulo(c(:) - 1, raxis(1:rdim)) + 1
162 this_corner = kl123(c(1), 1, 1)
163 this%simplices(n,ik) = this_corner
170 integer,
parameter :: submesh_triangles(2,3) = reshape([ &
172 1, 4, 3], shape(submesh_triangles), order=[2, 1])
174 integer,
parameter :: b(10,3) = reshape([ &
185 ], shape(b), order=[2, 1])
187 integer :: i, j, ip1, jp1, it, n
188 integer :: corners(4,2), v(3,2), c(2)
189 integer :: this_triangle(3), this_corner
191 this%n_simplices = 2 * npoints
193 this%sdim = merge(10, 3, opt)
194 safe_allocate(this%simplices(this%n_simplices, this%sdim))
198 ip1 = modulo(i, raxis(1)) + 1
199 jp1 = modulo(j, raxis(2)) + 1
200 corners(:,:) = reshape([ &
204 i , jp1 ], shape(corners), order=[2, 1])
206 do it = 1,
size(submesh_triangles, 1)
207 n = (it - 1) + 2 * ((j - 1) + raxis(2) * (i - 1)) + 1
208 this_triangle(:) = submesh_triangles(it, :)
209 v(1,:) = corners(this_triangle(1), :)
210 v(2,:) = corners(this_triangle(2), :)
211 v(3,:) = corners(this_triangle(3), :)
213 c(:) = b(ik,1) * v(1,:) + b(ik,2) * v(2,:) + b(ik,3) * v(3,:)
214 c(:) = modulo(c(:) - 1, raxis(1:rdim)) + 1
215 this_corner = kl123(c(1), c(2), 1)
216 this%simplices(n,ik) = this_corner
224 integer,
parameter :: submesh_tetras(6,4) = reshape([ &
230 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
232 integer,
parameter :: b(20,4) = reshape([ &
253 ], shape(b), order=[2, 1])
255 integer :: i, j, k, ip1, jp1, kp1, it, n
256 integer :: corners(8,3), v(4,3), c(3)
257 integer :: this_tetra(4), this_corner
259 this%n_simplices = 6 * npoints
261 this%sdim = merge(20, 4, opt)
262 safe_allocate(this%simplices(this%n_simplices, this%sdim))
267 ip1 = modulo(i, raxis(1)) + 1
268 jp1 = modulo(j, raxis(2)) + 1
269 kp1 = modulo(k, raxis(3)) + 1
270 corners(:,:) = reshape([ &
278 ip1 , jp1 , kp1 ], shape(corners), order=[2, 1])
280 do it = 1,
size(submesh_tetras, 1)
281 n = (it - 1) + 6 * ((k - 1) + raxis(3) * ((j - 1) + raxis(2) * (i - 1))) + 1
282 this_tetra(:) = submesh_tetras(it, :)
283 v(1,:) = corners(this_tetra(1), :)
284 v(2,:) = corners(this_tetra(2), :)
285 v(3,:) = corners(this_tetra(3), :)
286 v(4,:) = corners(this_tetra(4), :)
288 c(:) = b(ik,1) * v(1,:) + b(ik,2) * v(2,:) + b(ik,3) * v(3,:) + b(ik,4) * v(4,:)
289 c(:) = modulo(c(:) - 1, raxis(1:rdim)) + 1
290 this_corner = kl123(c(1), c(2), c(3))
291 this%simplices(n,ik) = this_corner
300 safe_deallocate_a(kl123)
311 safe_deallocate_a(this%simplices)
323 integer,
intent(in) :: rdim
324 real(real64),
intent(in) :: esimplex(:)
325 real(real64),
intent(in) :: ef
326 real(real64),
intent(out) :: weights(:)
327 real(real64),
intent(out) :: dos
350 integer,
intent(in) :: rdim
351 real(real64),
intent(in) :: esimplex(:)
352 real(real64),
intent(in) :: ef
353 real(real64),
intent(out) :: dos
376 real(real64),
intent(in) :: esegment(:)
377 real(real64),
intent(in) :: ef
378 real(real64),
intent(out) :: weights(:)
379 real(real64),
intent(out) :: dos
381 real(real64) :: e(2), e1, e2
385 real(real64),
parameter :: vt_vg = 1.0_real64
386 real(real64),
parameter :: vt_2vg = vt_vg / 2.0_real64
388 real(real64),
parameter :: p(2,4) = 1.0_real64 / 60.0_real64 * reshape([ &
390 1 , 64 , -2 , -3 ], shape(p), order=[2, 1])
394 select case (
size(esegment))
401 do i = 1,
size(esegment)
402 e(:) = e(:) + p(:,i) * esegment(i)
417 elseif (e2 < ef)
then
420 elseif (e1 < ef .and. ef <= e2)
then
422 real(real64) :: e21, c
424 c = vt_2vg * (ef - e1) / e21
427 2.0_real64 - (ef - e1) / e21, &
436 select case (
size(esegment))
452 real(real64),
intent(in) :: esegment(:)
453 real(real64),
intent(in) :: ef
454 real(real64),
intent(out) :: dos
456 real(real64) :: e(2), e1, e2
458 real(real64),
parameter :: vt_vg = 1.0_real64
460 real(real64),
parameter :: p(2,4) = 1.0_real64 / 60.0_real64 * reshape([ &
462 1 , 64 , -2 , -3 ], shape(p), order=[2, 1])
466 select case (
size(esegment))
473 do i = 1,
size(esegment)
474 e(:) = e(:) + p(:,i) * esegment(i)
485 if (ef <= e1 .or. e2 < ef)
then
487 elseif (e1 < ef .and. ef <= e2)
then
508 real(real64),
intent(in) :: etriangle(:)
509 real(real64),
intent(in) :: ef
510 real(real64),
intent(out) :: weights(:)
511 real(real64),
intent(out) :: dos
513 real(real64) :: e(3), e1, e2, e3
517 real(real64),
parameter :: vt_vg = 1.0_real64 / 2.0_real64
518 real(real64),
parameter :: vt_3vg = vt_vg / 3.0_real64
520 real(real64),
parameter :: p(3,10) = 1.0_real64 / 360.0_real64 * reshape([ &
521 402 , 0 , 6 , -13 , 5 , -17 , -13 , -11 , 7 , -6 , &
522 6 , 396 , 6 , -9 , -15 , 3 , 3 , -15 , -9 , -6 , &
523 6 , 0 , 402 , 7 , -11 , -13 , -17 , 5 , -13 , -6 &
524 ], shape(p), order=[2, 1])
528 select case (
size(etriangle))
535 do i = 1,
size(etriangle)
536 e(:) = e(:) + p(:,i) * etriangle(i)
552 elseif (e3 < ef)
then
555 elseif (e1 < ef .and. ef <= e2)
then
557 real(real64) :: e21, e31, c
560 c = vt_3vg * (ef - e1) ** 2 / (e21 * e31)
563 3.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31), &
567 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
569 elseif (e2 < ef .and. ef <= e3)
then
571 real(real64) :: e23, e31, c1, c2
575 c2 = vt_3vg * (ef - e3) ** 2 / (e23 * e31)
578 c1 - c2 * (ef - e3) / e31, &
579 c1 + c2 * (ef - e3) / e23, &
580 c1 + c2 * (3.0_real64 - (ef - e3) * (m_one / e23 - m_one / e31))]
582 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
588 select case (
size(etriangle))
590 weights(idx) = w + dos * (sum(e) - 3.0_real64 * e) / 24.0_real64
604 real(real64),
intent(in) :: etriangle(:)
605 real(real64),
intent(in) :: ef
606 real(real64),
intent(out) :: dos
608 real(real64) :: e(3), e1, e2, e3
610 real(real64),
parameter :: vt_vg = 1.0_real64 / 2.0_real64
612 real(real64),
parameter :: p(3,10) = 1.0_real64 / 360.0_real64 * reshape([ &
613 402 , 0 , 6 , -13 , 5 , -17 , -13 , -11 , 7 , -6 , &
614 6 , 396 , 6 , -9 , -15 , 3 , 3 , -15 , -9 , -6 , &
615 6 , 0 , 402 , 7 , -11 , -13 , -17 , 5 , -13 , -6 &
616 ], shape(p), order=[2, 1])
620 select case (
size(etriangle))
627 do i = 1,
size(etriangle)
628 e(:) = e(:) + p(:,i) * etriangle(i)
640 if (ef <= e1 .or. e3 < ef)
then
642 elseif (e1 < ef .and. ef <= e2)
then
644 real(real64) :: e21, e31
647 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
649 elseif (e2 < ef .and. ef <= e3)
then
651 real(real64) :: e23, e31
654 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
674 real(real64),
intent(in) :: etetra(:)
675 real(real64),
intent(in) :: ef
676 real(real64),
intent(out) :: weights(:)
677 real(real64),
intent(out) :: dos
679 real(real64) :: e(4), e1, e2, e3, e4
683 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
684 real(real64),
parameter :: vt_4vg = vt_vg / 4.0_real64
686 real(real64),
parameter :: p(4,20) = 1.0_real64 / 1260.0_real64 * reshape([ &
687 1440, 0 , 30 , 0 , -38 , 7 , 17 , -28 , -56 , 9 , -46 , 9 , -38 , -28 , 17 , 7 , -18 , -18 , 12 , -18 , &
688 0 , 1440, 0 , 30 , -28 , -38 , 7 , 17 , 9 , -56 , 9 , -46 , 7 , -38 , -28 , 17 , -18 , -18 , -18 , 12 , &
689 30 , 0 , 1440, 0 , 17 , -28 , -38 , 7 , -46 , 9 , -56 , 9 , 17 , 7 , -38 , -28 , 12 , -18 , -18 , -18 , &
690 0 , 30 , 0 , 1440, 7 , 17 , -28 , -38 , 9 , -46 , 9 , -56 , -28 , 17 , 7 , -38 , -18 , 12 , -18 , -18 &
691 ], shape(p), order=[2, 1])
695 select case (
size(etetra))
702 do i = 1,
size(etetra)
703 e(:) = e(:) + p(:,i) * etetra(i)
720 elseif (e4 < ef)
then
723 elseif (e1 < ef .and. ef <= e2)
then
725 real(real64) :: e21, e31, e41, c
729 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
732 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
737 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
739 elseif (e2 < ef .and. ef <= e3)
then
741 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
747 c1 = vt_4vg * (ef - e1) ** 2 / (e41 * e31)
748 c2 = vt_4vg * (ef - e1) * (ef - e2) * (e3 - ef) / (e41 * e32 * e31)
749 c3 = vt_4vg * (ef - e2) ** 2 * (e4 - ef) / (e42 * e32 * e41)
752 c1 + (c1 + c2) * (e3 - ef) / e31 + (c1 + c2 + c3) * (e4 - ef) / e41, &
753 c1 + c2 + c3 + (c2 + c3) * (e3 - ef) / e32 + c3 * (e4 - ef) / e42, &
754 (c1 + c2) * (ef - e1) / e31 + (c2 + c3) * (ef - e2) / e32, &
755 (c1 + c2 + c3) * (ef - e1) / e41 + c3 * (ef - e2) / e42]
757 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
758 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
760 elseif (e3 < ef .and. ef <= e4)
then
762 real(real64) :: e41, e42, e43, c
766 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
768 w(:) = vt_4vg - c * [ &
772 4.0_real64 - (e4 - ef) * (m_one / e41 + m_one / e42 + m_one / e43)]
774 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
780 select case (
size(etetra))
782 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
802 real(real64),
intent(in) :: etetra(:)
803 real(real64),
intent(in) :: ef
804 real(real64),
intent(out) :: dos
806 real(real64) :: e(4), e1, e2, e3, e4
808 real(real64),
parameter :: vt_vg = 1.0_real64 / 6.0_real64
810 real(real64),
parameter :: p(4,20) = 1.0_real64 / 1260.0_real64 * reshape([ &
811 1440, 0 , 30 , 0 , -38 , 7 , 17 , -28 , -56 , 9 , -46 , 9 , -38 , -28 , 17 , 7 , -18 , -18 , 12 , -18 , &
812 0 , 1440, 0 , 30 , -28 , -38 , 7 , 17 , 9 , -56 , 9 , -46 , 7 , -38 , -28 , 17 , -18 , -18 , -18 , 12 , &
813 30 , 0 , 1440, 0 , 17 , -28 , -38 , 7 , -46 , 9 , -56 , 9 , 17 , 7 , -38 , -28 , 12 , -18 , -18 , -18 , &
814 0 , 30 , 0 , 1440, 7 , 17 , -28 , -38 , 9 , -46 , 9 , -56 , -28 , 17 , 7 , -38 , -18 , 12 , -18 , -18 &
815 ], shape(p), order=[2, 1])
819 select case (
size(etetra))
826 do i = 1,
size(etetra)
827 e(:) = e(:) + p(:,i) * etetra(i)
840 if (e1 >= ef .or. e4 < ef)
then
842 elseif (e1 < ef .and. ef <= e2)
then
844 real(real64) :: e21, e31, e41
848 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
850 elseif (e2 < ef .and. ef <= e3)
then
852 real(real64) :: e21, e31, e32, e41, e42
858 dos = vt_vg / (e31 * e41) * (3.0_real64 * e21 + 6.0_real64 * (ef - e2) &
859 - 3.0_real64 * (e31 + e42) * (ef - e2) ** 2 / (e32 * e42))
861 elseif (e3 < ef .and. ef <= e4)
then
863 real(real64) :: e41, e42, e43
867 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(rdim, esimplex, eF, weights, dos)
Get the weights and DOS contribution of a single simplex.
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_dos_1d(esegment, eF, dos)
Get only the DOS contribution of a single segment.
subroutine, public simplex_end(this)
Destructor for linear simplex methods.
subroutine, public simplex_weights_1d(esegment, eF, weights, dos)
Get the weights and DOS contribution of a single segment.
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.
subroutine, public simplex_dos(rdim, esimplex, eF, dos)
Get only the DOS contribution of a single simplex.
This module is intended to contain "only mathematical" functions and procedures.