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, &
46
47 type simplex_t
48 integer :: dim
49 integer :: n_simplices
50 integer, allocatable :: simplices(:,:)
51 contains
52 final :: simplex_end
53 end type simplex_t
54
55 interface simplex_t
56 module procedure simplex_init
57 end interface simplex_t
58
59contains
60
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
83
84 real(real64) :: kmin(dim)
85 integer :: ik, npoints
86
87 integer :: ix(dim)
88 integer, allocatable :: kl123(:,:,:)
89
90 integer :: rdim, raxis(3)
91
92 push_sub(simplex_init)
93
94 if (nshifts /= 1) then
95 message(1) = "The linear tetrahedron method only works for automatic k-point grids with a single shift"
96 call messages_fatal(1)
97 end if
98
99 safe_allocate(this)
100 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
101
102 npoints = product(naxis)
103 kmin = minval(kpoints, 2)
104
105 do ik = 1, npoints
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)
109 end do
110
111 rdim = sum(merge(1, 0, naxis > 1))
112 raxis(1:rdim) = pack(naxis, naxis > 1)
113
114 if (any(raxis(1:rdim) /= naxis(1:rdim))) then
115 message(1) = "The periodic dimensions must be consecutive"
117 end if
118
119 select case (rdim)
120 case default
121 assert(.false.)
122 case (1)
123 call messages_not_implemented("The linear tetrahedron for 1D systems")
124 case (2)
125 block
126 integer, parameter :: submesh_triangles(2,3) = reshape([ &
127 1, 2, 3, &
128 1, 4, 3], shape(submesh_triangles), order=[2, 1])
129
130 integer :: i, j, ip1, jp1, it, n
131 integer :: corners(4,2)
132 integer :: this_triangle(3), this_corner
133
134 this%n_simplices = 2 * npoints
135 this%dim = 3
136 safe_allocate(this%simplices(this%n_simplices, this%dim))
137
138 do i = 1, raxis(1)
139 do j = 1, raxis(2)
140 ip1 = mod(i, raxis(1)) + 1
141 jp1 = mod(j, raxis(2)) + 1
142 corners(:,:) = reshape([ &
143 i , j , &
144 ip1 , j , &
145 ip1 , jp1 , &
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, :)
151 do ik = 1, 3
152 this_corner = kl123(corners(this_triangle(ik), 1), corners(this_triangle(ik), 2), 1)
153 this%simplices(n,ik) = this_corner
154 end do
155 end do
156 end do
157 end do
158 end block
159 case(3)
160 block
161 integer, parameter :: submesh_tetras(6,4) = reshape([ &
162 1, 2, 3, 6, &
163 1, 3, 5, 6, &
164 3, 5, 6, 7, &
165 3, 6, 7, 8, &
166 3, 4, 6, 8, &
167 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
168
169 integer :: i, j, k, ip1, jp1, kp1, it, n
170 integer :: corners(8,3)
171 integer :: this_tetra(4), this_corner
172
173 this%n_simplices = 6 * npoints
174 this%dim = 4
175 safe_allocate(this%simplices(this%n_simplices, this%dim))
176
177 do i = 1, raxis(1)
178 do j = 1, raxis(2)
179 do k = 1, raxis(3)
180 ip1 = mod(i, raxis(1)) + 1
181 jp1 = mod(j, raxis(2)) + 1
182 kp1 = mod(k, raxis(3)) + 1
183 corners(:,:) = reshape([ &
184 i , j , k , &
185 ip1 , j , k , &
186 i , jp1 , k , &
187 ip1 , jp1 , k , &
188 i , j , kp1 , &
189 ip1 , j , kp1 , &
190 i , jp1 , kp1 , &
191 ip1 , jp1 , kp1 ], shape(corners), order=[2, 1])
192
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, :)
196 do ik = 1, 4
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
199 end do
200 end do
201 end do
202 end do
203 end do
204 end block
205 end select
206
207 safe_deallocate_a(kl123)
208
209 pop_sub(simplex_init)
210 end function simplex_init
211
215 subroutine simplex_end(this)
216 type(simplex_t), intent(inout) :: this
217 push_sub(simplex_end)
218 safe_deallocate_a(this%simplices)
219 pop_sub(simplex_end)
220 end subroutine simplex_end
221
231 subroutine simplex_weights_2d(etriangle, eF, weights, dos)
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
236
237 real(real64) :: e(3), e1, e2, e3
238 real(real64) :: w(3)
239 integer :: idx(3)
240
241 real(real64), parameter :: vt_vg = 1.0_real64 / 2.0_real64
242 real(real64), parameter :: vt_3vg = vt_vg / 3.0_real64
243
244 ! no PUSH_SUB, called too often
245
246 idx = [1,2,3]
247 e = etriangle
248 call sort(e, idx)
249 e1 = e(1)
250 e2 = e(2)
251 e3 = e(3)
252
253 if (ef <= e1) then
254 w(:) = m_zero
255 dos = m_zero
256 elseif (e3 < ef) then
257 w(:) = vt_3vg
258 dos = m_zero
259 elseif (e1 < ef .and. ef <= e2) then
260 block
261 real(real64) :: e21, e31, c
262 e21 = e2 - e1
263 e31 = e3 - e1
264 c = vt_3vg * (ef - e1) ** 2 / (e21 * e31)
265
266 w(:) = c * [ &
267 3.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31), &
268 (ef - e1) / e21, &
269 (ef - e1) / e31]
270
271 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
272 end block
273 elseif (e2 < ef .and. ef <= e3) then
274 block
275 real(real64) :: e23, e31, c1, c2
276 e23 = e2 - e3
277 e31 = e3 - e1
278 c1 = vt_3vg
279 c2 = vt_3vg * (ef - e3) ** 2 / (e23 * e31)
280
281 w(:) = [ &
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))]
285
286 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
287 end block
288 else
289 assert(.false.)
290 end if
291
292 weights(idx) = w + dos * (sum(e) - 3.0_real64 * e) / 24.0_real64
293 end subroutine simplex_weights_2d
294
300 subroutine simplex_dos_2d(etriangle, eF, dos)
301 real(real64), intent(in) :: etriangle(3)
302 real(real64), intent(in) :: ef
303 real(real64), intent(out) :: dos
304
305 real(real64) :: e(3), e1, e2, e3
306
307 real(real64), parameter :: vt_vg = 1.0_real64 / 2.0_real64
308
309 ! no PUSH_SUB, called too often
311 e = etriangle
312 call sort(e)
313 e1 = e(1)
314 e2 = e(2)
315 e3 = e(3)
316
317 if (ef <= e1 .or. e3 < ef) then
318 dos = m_zero
319 elseif (e1 < ef .and. ef <= e2) then
320 block
321 real(real64) :: e21, e31
322 e21 = e2 - e1
323 e31 = e3 - e1
324 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
325 end block
326 elseif (e2 < ef .and. ef <= e3) then
327 block
328 real(real64) :: e23, e31
329 e23 = e2 - e3
330 e31 = e3 - e1
331 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
332 end block
333 else
334 assert(.false.)
335 end if
336 end subroutine simplex_dos_2d
337
347 subroutine simplex_weights_3d(etetra, eF, weights, dos)
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
352
353 real(real64) :: e(4), e1, e2, e3, e4
354 real(real64) :: w(4)
355 integer :: idx(4)
356
357 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
358 real(real64), parameter :: vt_4vg = vt_vg / 4.0_real64
359
360 ! no PUSH_SUB, called too often
361
362 idx = [1,2,3,4]
363 e = etetra
364 call sort(e, idx)
365 e1 = e(1)
366 e2 = e(2)
367 e3 = e(3)
368 e4 = e(4)
369
370 if (e1 >= ef) then
371 w(:) = m_zero
372 dos = m_zero
373 elseif (e4 < ef) then
374 w(:) = vt_4vg
375 dos = m_zero
376 elseif (e1 < ef .and. ef <= e2) then
377 block
378 real(real64) :: e21, e31, e41, c
379 e21 = e2 - e1
380 e31 = e3 - e1
381 e41 = e4 - e1
382 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
383
384 w(:) = c * [ &
385 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
386 (ef - e1) / e21, &
387 (ef - e1) / e31, &
388 (ef - e1) / e41]
389
390 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
391 end block
392 elseif (e2 < ef .and. ef <= e3) then
393 block
394 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
395 e21 = e2 - e1
396 e31 = e3 - e1
397 e32 = e3 - e2
398 e41 = e4 - e1
399 e42 = e4 - e2
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)
403
404 w(:) = [ &
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]
409
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))
412 end block
413 elseif (e3 < ef .and. ef <= e4) then
414 block
415 real(real64) :: e41, e42, e43, c
416 e41 = e4 - e1
417 e42 = e4 - e2
418 e43 = e4 - e3
419 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
420
421 w(:) = vt_4vg - c * [ &
422 (e4 - ef) / e41, &
423 (e4 - ef) / e42, &
424 (e4 - ef) / e43, &
425 4.0_real64 - (e4 - ef) * (m_one / e41 + m_one / e42 + m_one / e43)]
426
427 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
428 end block
429 else
430 assert(.false.)
431 end if
432
433 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
434 end subroutine simplex_weights_3d
435
444 subroutine simplex_dos_3d(etetra, eF, dos)
445 real(real64), intent(in) :: etetra(4)
446 real(real64), intent(in) :: ef
447 real(real64), intent(out) :: dos
448
449 real(real64) :: e(4), e1, e2, e3, e4
450
451 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
452
453 ! no PUSH_SUB, called too often
454
455 e = etetra
456 call sort(e)
457 e1 = e(1)
458 e2 = e(2)
459 e3 = e(3)
460 e4 = e(4)
461
462 if (e1 >= ef .or. e4 < ef) then
463 dos = m_zero
464 elseif (e1 < ef .and. ef <= e2) then
465 block
466 real(real64) :: e21, e31, e41
467 e21 = e2 - e1
468 e31 = e3 - e1
469 e41 = e4 - e1
470 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
471 end block
472 elseif (e2 < ef .and. ef <= e3) then
473 block
474 real(real64) :: e21, e31, e32, e41, e42
475 e21 = e2 - e1
476 e31 = e3 - e1
477 e32 = e3 - e2
478 e41 = e4 - e1
479 e42 = e4 - e2
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))
482 end block
483 elseif (e3 < ef .and. ef <= e4) then
484 block
485 real(real64) :: e41, e42, e43
486 e41 = e4 - e1
487 e42 = e4 - e2
488 e43 = e4 - e3
489 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
490 end block
491 else
492 assert(.false.)
493 end if
494 end subroutine simplex_dos_3d
495
496end module simplex_oct_m
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:156
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 messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
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_weights_2d(etriangle, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
Definition: simplex.F90:327
subroutine, public simplex_dos_3d(etetra, eF, dos)
Get only the DOS contribution of a single tetrahedron.
Definition: simplex.F90:540
subroutine, public simplex_end(this)
Destructor for linear simplex methods.
Definition: simplex.F90:311
type(simplex_t) function, pointer, public simplex_init(dim, naxis, nshifts, shift, kpoints, equiv)
Constructor for linear simplex methods.
Definition: simplex.F90:171
subroutine, public simplex_dos_2d(etriangle, eF, dos)
Get only the DOS contribution of a single triangle.
Definition: simplex.F90:396
subroutine, public simplex_weights_3d(etetra, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
Definition: simplex.F90:443
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119