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, &
50
51 type simplex_t
52 integer :: rdim
53 integer :: sdim
54 integer :: n_simplices
55 integer, allocatable :: simplices(:,:)
56 contains
57 final :: simplex_end
58 end type simplex_t
59
60 interface simplex_t
61 module procedure simplex_init
62 end interface simplex_t
63
64contains
65
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
90
91 real(real64) :: kmin(dim)
92 integer :: ik, npoints
93
94 integer :: ix(dim)
95 integer, allocatable :: kl123(:,:,:)
96
97 integer :: rdim, raxis(3)
98
99 push_sub(simplex_init)
100
101 if (nshifts /= 1) then
102 message(1) = "The linear tetrahedron method only works for automatic k-point grids with a single shift"
103 call messages_fatal(1)
104 end if
105
106 safe_allocate(this)
107 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
108
109 npoints = product(naxis)
110 kmin = minval(kpoints, 2)
111
112 do ik = 1, npoints
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)
116 end do
117
118 rdim = sum(merge(1, 0, naxis > 1))
119 raxis(1:rdim) = pack(naxis, naxis > 1)
120
121 if (any(raxis(1:rdim) /= naxis(1:rdim))) then
122 message(1) = "The periodic dimensions must be consecutive"
123 call messages_fatal(1)
124 end if
125
126 select case (rdim)
127 case default
128 assert(.false.)
129 case (1)
130 block
131 integer, parameter :: submesh_segments(1,2) = reshape([ &
132 1, 2 ], shape(submesh_segments), order=[2, 1])
133 ! coordinates of corners and neighbors in barycentric coordinates
134 integer, parameter :: b(4,2) = reshape([ &
135 1 , 0 , & ! k1
136 0 , 1 , & ! k2
137 2 , -1, & ! 2 * k1 - k2
138 -1, 2 & ! 2 * k2 - k1
139 ], shape(b), order=[2, 1])
140
141 integer :: i, ip1, it, n
142 integer :: corners(2,1), v(2,1), c(1)
143 integer :: this_segment(2), this_corner
144
145 this%n_simplices = npoints
146 this%rdim = rdim
147 this%sdim = merge(4, 2, opt)
148 safe_allocate(this%simplices(this%n_simplices, this%sdim))
150 do i = 1, raxis(1)
151 ip1 = modulo(i, raxis(1)) + 1
152 corners(:,:) = reshape([ i , ip1 ], shape(corners), order=[2, 1])
153
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), :)
159 do ik = 1, this%sdim
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
164 end do
165 end do
166 end do
167 end block
168 case (2)
169 block
170 integer, parameter :: submesh_triangles(2,3) = reshape([ &
171 1, 2, 3, &
172 1, 4, 3], shape(submesh_triangles), order=[2, 1])
173 ! coordinates of corners and neighbors in barycentric coordinates
174 integer, parameter :: b(10,3) = reshape([ &
175 1 , 0 , 0 , & ! k1
176 0 , 1 , 0 , & ! k2
177 0 , 0 , 1 , & ! k3
178 2 , -1, 0 , & ! 2 * k1 - k2
179 0 , 2 , -1, & ! 2 * k2 - k3
180 2 , 0 , -1, & ! 2 * k1 - k3
181 -1, 0 , 2 , & ! 2 * k3 - k1
182 -1, 2 , 0 , & ! 2 * k2 - k1
183 0 , -1, 2 , & ! 2 * k3 - k2
184 1 , -1, 1 & ! k1 - k2 + k3
185 ], shape(b), order=[2, 1])
186
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
190
191 this%n_simplices = 2 * npoints
192 this%rdim = rdim
193 this%sdim = merge(10, 3, opt)
194 safe_allocate(this%simplices(this%n_simplices, this%sdim))
195
196 do i = 1, raxis(1)
197 do j = 1, raxis(2)
198 ip1 = modulo(i, raxis(1)) + 1
199 jp1 = modulo(j, raxis(2)) + 1
200 corners(:,:) = reshape([ &
201 i , j , &
202 ip1 , j , &
203 ip1 , jp1 , &
204 i , jp1 ], shape(corners), order=[2, 1])
205
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), :)
212 do ik = 1, this%sdim
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
217 end do
218 end do
219 end do
220 end do
221 end block
222 case(3)
223 block
224 integer, parameter :: submesh_tetras(6,4) = reshape([ &
225 1, 2, 3, 6, &
226 1, 3, 5, 6, &
227 3, 5, 6, 7, &
228 3, 6, 7, 8, &
229 3, 4, 6, 8, &
230 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
231 ! coordinates of corners and neighbors in barycentric coordinates
232 integer, parameter :: b(20,4) = reshape([ &
233 1 , 0 , 0 , 0 , & ! k1
234 0 , 1 , 0 , 0 , & ! k2
235 0 , 0 , 1 , 0 , & ! k3
236 0 , 0 , 0 , 1 , & ! k4
237 2 , -1, 0 , 0 , & ! 2 * k1 - k2
238 0 , 2 , -1, 0 , & ! 2 * k2 - k3
239 0 , 0 , 2 , -1, & ! 2 * k3 - k4
240 -1, 0 , 0 , 2 , & ! 2 * k4 - k1
241 2 , 0 , -1, 0 , & ! 2 * k1 - k3
242 0 , 2 , 0 , -1, & ! 2 * k2 - k4
243 -1, 0 , 2 , 0 , & ! 2 * k3 - k1
244 0 , -1, 0 , 2 , & ! 2 * k4 - k2
245 2 , 0 , 0 , -1, & ! 2 * k1 - k4
246 -1, 2 , 0 , 0 , & ! 2 * k2 - k1
247 0 , -1, 2 , 0 , & ! 2 * k3 - k2
248 0 , 0 , -1, 2 , & ! 2 * k4 - k3
249 -1, 1 , 0 , 1 , & ! k4 - k1 + k2
250 1 , -1, 1 , 0 , & ! k1 - k2 + k3
251 0 , 1 , -1, 1 , & ! k2 - k3 + k4
252 1 , 0 , 1 , -1 & ! k3 - k4 + k1
253 ], shape(b), order=[2, 1])
254
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
258
259 this%n_simplices = 6 * npoints
260 this%rdim = rdim
261 this%sdim = merge(20, 4, opt)
262 safe_allocate(this%simplices(this%n_simplices, this%sdim))
263
264 do i = 1, raxis(1)
265 do j = 1, raxis(2)
266 do k = 1, raxis(3)
267 ip1 = modulo(i, raxis(1)) + 1
268 jp1 = modulo(j, raxis(2)) + 1
269 kp1 = modulo(k, raxis(3)) + 1
270 corners(:,:) = reshape([ &
271 i , j , k , &
272 ip1 , j , k , &
273 i , jp1 , k , &
274 ip1 , jp1 , k , &
275 i , j , kp1 , &
276 ip1 , j , kp1 , &
277 i , jp1 , kp1 , &
278 ip1 , jp1 , kp1 ], shape(corners), order=[2, 1])
279
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), :)
287 do ik = 1, this%sdim
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
292 end do
293 end do
294 end do
295 end do
296 end do
297 end block
298 end select
299
300 safe_deallocate_a(kl123)
301
302 pop_sub(simplex_init)
303 end function simplex_init
304
308 subroutine simplex_end(this)
309 type(simplex_t), intent(inout) :: this
310 push_sub(simplex_end)
311 safe_deallocate_a(this%simplices)
312 pop_sub(simplex_end)
313 end subroutine simplex_end
314
322 subroutine simplex_weights(rdim, esimplex, eF, weights, dos)
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
328
329 ! no PUSH_SUB, called too often
330
331 select case (rdim)
332 case (1)
333 call simplex_weights_1d(esimplex, ef, weights, dos)
334 case (2)
335 call simplex_weights_2d(esimplex, ef, weights, dos)
336 case (3)
337 call simplex_weights_3d(esimplex, ef, weights, dos)
338 case default
339 assert(.false.)
340 end select
341 end subroutine simplex_weights
342
349 subroutine simplex_dos(rdim, esimplex, eF, dos)
350 integer, intent(in) :: rdim
351 real(real64), intent(in) :: esimplex(:)
352 real(real64), intent(in) :: ef
353 real(real64), intent(out) :: dos
354
355 ! no PUSH_SUB, called too often
356
357 select case (rdim)
358 case (1)
359 call simplex_dos_1d(esimplex, ef, dos)
360 case (2)
361 call simplex_dos_2d(esimplex, ef, dos)
362 case (3)
363 call simplex_dos_3d(esimplex, ef, dos)
364 case default
365 assert(.false.)
366 end select
367 end subroutine simplex_dos
368
375 subroutine simplex_weights_1d(esegment, eF, weights, 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
380
381 real(real64) :: e(2), e1, e2
382 real(real64) :: w(2)
383 integer :: idx(2)
384
385 real(real64), parameter :: vt_vg = 1.0_real64
386 real(real64), parameter :: vt_2vg = vt_vg / 2.0_real64
387
388 real(real64), parameter :: p(2,4) = 1.0_real64 / 60.0_real64 * reshape([ &
389 64 , 1 , -3 , -2 , &
390 1 , 64 , -2 , -3 ], shape(p), order=[2, 1])
391
392 ! no PUSH_SUB, called too often
393
394 select case (size(esegment))
395 case (2)
396 e(:) = esegment(:)
397 case (4)
398 e(:) = m_zero
399 block
400 integer :: i
401 do i = 1, size(esegment)
402 e(:) = e(:) + p(:,i) * esegment(i)
403 end do
404 end block
405 case default
406 assert(.false.)
407 end select
408
409 idx = [1,2]
410 call sort(e, idx)
411 e1 = e(1)
412 e2 = e(2)
413
414 if (ef <= e1) then
415 w(:) = m_zero
416 dos = m_zero
417 elseif (e2 < ef) then
418 w(:) = vt_2vg
419 dos = m_zero
420 elseif (e1 < ef .and. ef <= e2) then
421 block
422 real(real64) :: e21, c
423 e21 = e2 - e1
424 c = vt_2vg * (ef - e1) / e21
425
426 w(:) = c * [ &
427 2.0_real64 - (ef - e1) / e21, &
428 (ef - e1) / e21]
429
430 dos = vt_vg / e21
431 end block
432 else
433 assert(.false.)
434 end if
435
436 select case (size(esegment))
437 case (2)
438 weights(idx) = w ! TODO: Bloechl correction
439 case (4)
440 weights(idx) = w ! no Bloechl correction
441 case default
442 assert(.false.)
443 end select
444 end subroutine simplex_weights_1d
445
451 subroutine simplex_dos_1d(esegment, eF, dos)
452 real(real64), intent(in) :: esegment(:)
453 real(real64), intent(in) :: ef
454 real(real64), intent(out) :: dos
455
456 real(real64) :: e(2), e1, e2
457
458 real(real64), parameter :: vt_vg = 1.0_real64
459
460 real(real64), parameter :: p(2,4) = 1.0_real64 / 60.0_real64 * reshape([ &
461 64 , 1 , -3 , -2 , &
462 1 , 64 , -2 , -3 ], shape(p), order=[2, 1])
463
464 ! no PUSH_SUB, called too often
465
466 select case (size(esegment))
467 case (2)
468 e(:) = esegment(:)
469 case (4)
470 e(:) = m_zero
471 block
472 integer :: i
473 do i = 1, size(esegment)
474 e(:) = e(:) + p(:,i) * esegment(i)
475 end do
476 end block
477 case default
478 assert(.false.)
479 end select
480
481 call sort(e)
482 e1 = e(1)
483 e2 = e(2)
484
485 if (ef <= e1 .or. e2 < ef) then
486 dos = m_zero
487 elseif (e1 < ef .and. ef <= e2) then
488 block
489 real(real64) :: e21
490 e21 = e2 - e1
491 dos = vt_vg / e21
492 end block
493 else
494 assert(.false.)
495 end if
496 end subroutine simplex_dos_1d
497
507 subroutine simplex_weights_2d(etriangle, eF, weights, dos)
508 real(real64), intent(in) :: etriangle(:)
509 real(real64), intent(in) :: ef
510 real(real64), intent(out) :: weights(:)
511 real(real64), intent(out) :: dos
512
513 real(real64) :: e(3), e1, e2, e3
514 real(real64) :: w(3)
515 integer :: idx(3)
516
517 real(real64), parameter :: vt_vg = 1.0_real64 / 2.0_real64
518 real(real64), parameter :: vt_3vg = vt_vg / 3.0_real64
519
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])
525
526 ! no PUSH_SUB, called too often
527
528 select case (size(etriangle))
529 case (3)
530 e(:) = etriangle(:)
531 case (10)
532 e(:) = m_zero
533 block
534 integer :: i
535 do i = 1, size(etriangle)
536 e(:) = e(:) + p(:,i) * etriangle(i)
537 end do
538 end block
539 case default
540 assert(.false.)
541 end select
542
543 idx = [1,2,3]
544 call sort(e, idx)
545 e1 = e(1)
546 e2 = e(2)
547 e3 = e(3)
548
549 if (ef <= e1) then
550 w(:) = m_zero
551 dos = m_zero
552 elseif (e3 < ef) then
553 w(:) = vt_3vg
554 dos = m_zero
555 elseif (e1 < ef .and. ef <= e2) then
556 block
557 real(real64) :: e21, e31, c
558 e21 = e2 - e1
559 e31 = e3 - e1
560 c = vt_3vg * (ef - e1) ** 2 / (e21 * e31)
561
562 w(:) = c * [ &
563 3.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31), &
564 (ef - e1) / e21, &
565 (ef - e1) / e31]
566
567 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
568 end block
569 elseif (e2 < ef .and. ef <= e3) then
570 block
571 real(real64) :: e23, e31, c1, c2
572 e23 = e2 - e3
573 e31 = e3 - e1
574 c1 = vt_3vg
575 c2 = vt_3vg * (ef - e3) ** 2 / (e23 * e31)
576
577 w(:) = [ &
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))]
581
582 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
583 end block
584 else
585 assert(.false.)
586 end if
587
588 select case (size(etriangle))
589 case (3)
590 weights(idx) = w + dos * (sum(e) - 3.0_real64 * e) / 24.0_real64
591 case (10)
592 weights(idx) = w ! no Bloechl correction
593 case default
594 assert(.false.)
595 end select
596 end subroutine simplex_weights_2d
597
603 subroutine simplex_dos_2d(etriangle, eF, dos)
604 real(real64), intent(in) :: etriangle(:)
605 real(real64), intent(in) :: ef
606 real(real64), intent(out) :: dos
607
608 real(real64) :: e(3), e1, e2, e3
609
610 real(real64), parameter :: vt_vg = 1.0_real64 / 2.0_real64
611
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])
617
618 ! no PUSH_SUB, called too often
619
620 select case (size(etriangle))
621 case (3)
622 e(:) = etriangle(:)
623 case (10)
624 e(:) = m_zero
625 block
626 integer :: i
627 do i = 1, size(etriangle)
628 e(:) = e(:) + p(:,i) * etriangle(i)
629 end do
630 end block
631 case default
632 assert(.false.)
633 end select
634
635 call sort(e)
636 e1 = e(1)
637 e2 = e(2)
638 e3 = e(3)
639
640 if (ef <= e1 .or. e3 < ef) then
641 dos = m_zero
642 elseif (e1 < ef .and. ef <= e2) then
643 block
644 real(real64) :: e21, e31
645 e21 = e2 - e1
646 e31 = e3 - e1
647 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
648 end block
649 elseif (e2 < ef .and. ef <= e3) then
650 block
651 real(real64) :: e23, e31
652 e23 = e2 - e3
653 e31 = e3 - e1
654 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
655 end block
656 else
657 assert(.false.)
658 end if
659 end subroutine simplex_dos_2d
660
673 subroutine simplex_weights_3d(etetra, eF, weights, dos)
674 real(real64), intent(in) :: etetra(:)
675 real(real64), intent(in) :: ef
676 real(real64), intent(out) :: weights(:)
677 real(real64), intent(out) :: dos
678
679 real(real64) :: e(4), e1, e2, e3, e4
680 real(real64) :: w(4)
681 integer :: idx(4)
682
683 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
684 real(real64), parameter :: vt_4vg = vt_vg / 4.0_real64
685
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])
692
693 ! no PUSH_SUB, called too often
694
695 select case (size(etetra))
696 case (4)
697 e(:) = etetra(:)
698 case (20)
699 e(:) = m_zero
700 block
701 integer :: i
702 do i = 1, size(etetra)
703 e(:) = e(:) + p(:,i) * etetra(i)
704 end do
705 end block
706 case default
707 assert(.false.)
708 end select
709
710 idx = [1,2,3,4]
711 call sort(e, idx)
712 e1 = e(1)
713 e2 = e(2)
714 e3 = e(3)
715 e4 = e(4)
716
717 if (e1 >= ef) then
718 w(:) = m_zero
719 dos = m_zero
720 elseif (e4 < ef) then
721 w(:) = vt_4vg
722 dos = m_zero
723 elseif (e1 < ef .and. ef <= e2) then
724 block
725 real(real64) :: e21, e31, e41, c
726 e21 = e2 - e1
727 e31 = e3 - e1
728 e41 = e4 - e1
729 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
730
731 w(:) = c * [ &
732 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
733 (ef - e1) / e21, &
734 (ef - e1) / e31, &
735 (ef - e1) / e41]
736
737 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
738 end block
739 elseif (e2 < ef .and. ef <= e3) then
740 block
741 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
742 e21 = e2 - e1
743 e31 = e3 - e1
744 e32 = e3 - e2
745 e41 = e4 - e1
746 e42 = e4 - e2
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)
750
751 w(:) = [ &
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]
756
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))
759 end block
760 elseif (e3 < ef .and. ef <= e4) then
761 block
762 real(real64) :: e41, e42, e43, c
763 e41 = e4 - e1
764 e42 = e4 - e2
765 e43 = e4 - e3
766 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
767
768 w(:) = vt_4vg - c * [ &
769 (e4 - ef) / e41, &
770 (e4 - ef) / e42, &
771 (e4 - ef) / e43, &
772 4.0_real64 - (e4 - ef) * (m_one / e41 + m_one / e42 + m_one / e43)]
773
774 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
775 end block
776 else
777 assert(.false.)
778 end if
779
780 select case (size(etetra))
781 case (4)
782 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
783 case (20)
784 weights(idx) = w ! no Bloechl correction
785 case default
786 assert(.false.)
787 end select
788 end subroutine simplex_weights_3d
789
801 subroutine simplex_dos_3d(etetra, eF, dos)
802 real(real64), intent(in) :: etetra(:)
803 real(real64), intent(in) :: ef
804 real(real64), intent(out) :: dos
805
806 real(real64) :: e(4), e1, e2, e3, e4
807
808 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
809
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])
816
817 ! no PUSH_SUB, called too often
818
819 select case (size(etetra))
820 case (4)
821 e(:) = etetra(:)
822 case (20)
823 e(:) = m_zero
824 block
825 integer :: i
826 do i = 1, size(etetra)
827 e(:) = e(:) + p(:,i) * etetra(i)
828 end do
829 end block
830 case default
831 assert(.false.)
832 end select
833
834 call sort(e)
835 e1 = e(1)
836 e2 = e(2)
837 e3 = e(3)
838 e4 = e(4)
839
840 if (e1 >= ef .or. e4 < ef) then
841 dos = m_zero
842 elseif (e1 < ef .and. ef <= e2) then
843 block
844 real(real64) :: e21, e31, e41
845 e21 = e2 - e1
846 e31 = e3 - e1
847 e41 = e4 - e1
848 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
849 end block
850 elseif (e2 < ef .and. ef <= e3) then
851 block
852 real(real64) :: e21, e31, e32, e41, e42
853 e21 = e2 - e1
854 e31 = e3 - e1
855 e32 = e3 - e2
856 e41 = e4 - e1
857 e42 = e4 - e2
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))
860 end block
861 elseif (e3 < ef .and. ef <= e4) then
862 block
863 real(real64) :: e41, e42, e43
864 e41 = e4 - e1
865 e42 = e4 - e2
866 e43 = e4 - e3
867 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
868 end block
869 else
870 assert(.false.)
871 end if
872 end subroutine simplex_dos_3d
873
874end 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:514
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:264
logical pure function, public not_in_openmp()
Definition: global.F90:530
character(len=100), public global_alloc_errmsg
Definition: global.F90:265
integer, public global_alloc_err
Definition: global.F90:263
real(real64), parameter, public m_one
Definition: global.F90:192
subroutine, public alloc_error(size, file, line)
Definition: messages.F90:668
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public dealloc_error(size, file, line)
Definition: messages.F90:679
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
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(rdim, esimplex, eF, weights, dos)
Get the weights and DOS contribution of a single simplex.
Definition: simplex.F90:418
subroutine, public simplex_weights_2d(etriangle, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
Definition: simplex.F90:603
subroutine, public simplex_dos_3d(etetra, eF, dos)
Get only the DOS contribution of a single tetrahedron.
Definition: simplex.F90:897
type(simplex_t) function, pointer, public simplex_init(dim, naxis, nshifts, shift, kpoints, equiv, opt)
Constructor for linear simplex methods.
Definition: simplex.F90:177
subroutine, public simplex_dos_1d(esegment, eF, dos)
Get only the DOS contribution of a single segment.
Definition: simplex.F90:547
subroutine, public simplex_end(this)
Destructor for linear simplex methods.
Definition: simplex.F90:404
subroutine, public simplex_weights_1d(esegment, eF, weights, dos)
Get the weights and DOS contribution of a single segment.
Definition: simplex.F90:471
subroutine, public simplex_dos_2d(etriangle, eF, dos)
Get only the DOS contribution of a single triangle.
Definition: simplex.F90:699
subroutine, public simplex_weights_3d(etetra, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
Definition: simplex.F90:769
subroutine, public simplex_dos(rdim, esimplex, eF, dos)
Get only the DOS contribution of a single simplex.
Definition: simplex.F90:445
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119