Octopus
simplex.F90
Go to the documentation of this file.
1!! Copyright (C) 2025 Octopus Developers
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module simplex_oct_m
22 use, intrinsic :: iso_fortran_env, only: real64
23 use debug_oct_m, only: debug
24#ifndef NDEBUG
26#endif
27 use global_oct_m, only: &
30 use profiling_oct_m, only: &
32 use sort_oct_m, only: sort
33
34 implicit none
35
36 private
37
38 public :: &
39 simplex_t, &
44
45 type simplex_t
46 integer :: dim
47 integer :: n_simplices
48 integer, allocatable :: simplices(:,:)
49 contains
50 final :: simplex_end
51 end type simplex_t
52
53 interface simplex_t
54 module procedure simplex_init
55 end interface simplex_t
56
57 integer, parameter :: submesh_tetras(6,4) = reshape([ &
58 1, 2, 3, 6, &
59 1, 3, 5, 6, &
60 3, 5, 6, 7, &
61 3, 6, 7, 8, &
62 3, 4, 6, 8, &
63 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
64
65contains
66
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
89
90 real(real64) :: kmin(dim)
91 integer :: ik, npoints
92
93 integer :: ix(dim)
94 integer, allocatable :: kl123(:,:,:)
95
96 push_sub(simplex_init)
97
98 if (dim /= 3 .or. .not. all(naxis > 1)) then
99 message(1) = "The linear tetrahedron method is currently only implemented for 3D systems"
100 call messages_fatal(1)
101 end if
102 if (nshifts /= 1) then
103 message(1) = "The linear tetrahedron method only works for automatic k-point grids with a single shift"
104 call messages_fatal(1)
105 end if
106
107 safe_allocate(this)
108 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
109
110 npoints = product(naxis)
111 kmin = minval(kpoints, 2)
112
113 do ik = 1, npoints
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)
117 end do
118
119 this%n_simplices = 6 * npoints
120 this%dim = dim
121 safe_allocate(this%simplices(this%n_simplices, 4))
122
123 block
124 integer :: i, j, k, ip1, jp1, kp1, itetra, n
125 integer :: corners(8,3)
126 integer :: this_tetra(4), this_corner
127 do i = 1, naxis(1)
128 do j = 1, naxis(2)
129 do k = 1, naxis(3)
130 ip1 = mod(i, naxis(1)) + 1
131 jp1 = mod(j, naxis(2)) + 1
132 kp1 = mod(k, naxis(3)) + 1
133 corners(:,:) = reshape([ &
134 i , j , k , &
135 ip1 , j , k , &
136 i , jp1 , k , &
137 ip1 , jp1 , k , &
138 i , j , kp1 , &
139 ip1 , j , kp1 , &
140 i , jp1 , kp1 , &
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, :)
146 do ik = 1, 4
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
149 end do
150 end do
151 end do
152 end do
153 end do
154 end block
155
156 safe_deallocate_a(kl123)
157
158 pop_sub(simplex_init)
159 end function simplex_init
160
164 subroutine simplex_end(this)
165 type(simplex_t), intent(inout) :: this
166 push_sub(simplex_end)
167 safe_deallocate_a(this%simplices)
168 pop_sub(simplex_end)
169 end subroutine simplex_end
170
180 subroutine simplex_weights_3d(etetra, eF, weights, dos)
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
185
186 real(real64) :: e(4), e1, e2, e3, e4
187 real(real64) :: w(4)
188 integer :: idx(4)
189
190 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
191 real(real64), parameter :: vt_4vg = vt_vg / 4.0_real64
192
193 ! no PUSH_SUB, called too often
194
195 idx = [1,2,3,4]
196 e = etetra
197 call sort(e, idx)
198 e1 = e(1)
199 e2 = e(2)
200 e3 = e(3)
201 e4 = e(4)
202
203 if (e1 >= ef) then
204 w(:) = m_zero
205 dos = m_zero
206 elseif (e4 < ef) then
207 w(:) = vt_4vg
208 dos = m_zero
209 elseif (e1 < ef .and. ef <= e2) then
210 block
211 real(real64) :: e21, e31, e41, c
212 e21 = e2 - e1
213 e31 = e3 - e1
214 e41 = e4 - e1
215 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
216
217 w(:) = c * [ &
218 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
219 (ef - e1) / e21, &
220 (ef - e1) / e31, &
221 (ef - e1) / e41]
222
223 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
224 end block
225 elseif (e2 < ef .and. ef <= e3) then
226 block
227 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
228 e21 = e2 - e1
229 e31 = e3 - e1
230 e32 = e3 - e2
231 e41 = e4 - e1
232 e42 = e4 - e2
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)
236
237 w(:) = [ &
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]
242
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))
245 end block
246 elseif (e3 < ef .and. ef <= e4) then
247 block
248 real(real64) :: e41, e42, e43, c
249 e41 = e4 - e1
250 e42 = e4 - e2
251 e43 = e4 - e3
252 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
253
254 w(:) = vt_4vg - c * [ &
255 (e4 - ef) / e41, &
256 (e4 - ef) / e42, &
257 (e4 - ef) / e43, &
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)
261 end block
262 else
263 assert(.false.)
264 end if
265
266 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
267 end subroutine simplex_weights_3d
268
277 subroutine simplex_dos_3d(etetra, eF, dos)
278 real(real64), intent(in) :: etetra(4)
279 real(real64), intent(in) :: ef
280 real(real64), intent(out) :: dos
281
282 real(real64) :: e(4), e1, e2, e3, e4
283
284 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
285
286 ! no PUSH_SUB, called too often
287
288 e = etetra
289 call sort(e)
290 e1 = e(1)
291 e2 = e(2)
292 e3 = e(3)
293 e4 = e(4)
294
295 if (e1 >= ef .or. e4 < ef) then
296 dos = m_zero
297 elseif (e1 < ef .and. ef <= e2) then
298 block
299 real(real64) :: e21, e31, e41
300 e21 = e2 - e1
301 e31 = e3 - e1
302 e41 = e4 - e1
303 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
304 end block
305 elseif (e2 < ef .and. ef <= e3) then
306 block
307 real(real64) :: e21, e31, e32, e41, e42
308 e21 = e2 - e1
309 e31 = e3 - e1
310 e32 = e3 - e2
311 e41 = e4 - e1
312 e42 = e4 - e2
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))
315 end block
316 elseif (e3 < ef .and. ef <= e4) then
317 block
318 real(real64) :: e41, e42, e43
319 e41 = e4 - e1
320 e42 = e4 - e2
321 e43 = e4 - e3
322 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
323 end block
324 else
325 assert(.false.)
326 end if
327 end subroutine simplex_dos_3d
328
329end module simplex_oct_m
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:151
type(debug_t), save, public debug
Definition: debug.F90:158
subroutine, public debug_pop_sub(sub_name)
Pop a routine from the debug trace.
Definition: debug.F90:513
subroutine, public debug_push_sub(sub_name)
Push a routine to the debug trace.
Definition: debug.F90:441
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
integer(int64), public global_sizeof
Definition: global.F90:260
logical pure function, public not_in_openmp()
Definition: global.F90:494
character(len=100), public global_alloc_errmsg
Definition: global.F90:261
integer, public global_alloc_err
Definition: global.F90:259
real(real64), parameter, public m_one
Definition: global.F90:192
subroutine, public alloc_error(size, file, line)
Definition: messages.F90:669
subroutine, public dealloc_error(size, file, line)
Definition: messages.F90:680
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
type(profile_vars_t), target, save, public prof_vars
Definition: profiling.F90:248
integer, parameter, public profiling_memory
Definition: profiling.F90:208
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public profiling_memory_deallocate(var, file, line, size)
Definition: profiling.F90:1404
subroutine, public profiling_memory_allocate(var, file, line, size_)
Definition: profiling.F90:1333
subroutine, public simplex_dos_3d(etetra, eF, dos)
Get only the DOS contribution of a single tetrahedron.
Definition: simplex.F90:373
subroutine, public simplex_end(this)
Destructor for linear simplex methods.
Definition: simplex.F90:260
type(simplex_t) function, pointer, public simplex_init(dim, naxis, nshifts, shift, kpoints, equiv)
Constructor for linear simplex methods.
Definition: simplex.F90:177
subroutine, public simplex_weights_3d(etetra, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
Definition: simplex.F90:276
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119