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
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
85
86 real(real64) :: kmin(dim)
87 integer :: ik, npoints
88
89 integer :: ix(dim)
90 integer, allocatable :: kl123(:,:,:)
91
92 integer :: rdim, raxis(3)
93
94 push_sub(simplex_init)
95
96 if (nshifts /= 1) then
97 message(1) = "The linear tetrahedron method only works for automatic k-point grids with a single shift"
98 call messages_fatal(1)
99 end if
100
101 safe_allocate(this)
102 safe_allocate_source(kl123(1:naxis(1), 1:naxis(2), 1:naxis(3)), -1)
103
104 npoints = product(naxis)
105 kmin = minval(kpoints, 2)
106
107 do ik = 1, npoints
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)
111 end do
112
113 rdim = sum(merge(1, 0, naxis > 1))
114 raxis(1:rdim) = pack(naxis, naxis > 1)
115
116 if (any(raxis(1:rdim) /= naxis(1:rdim))) then
117 message(1) = "The periodic dimensions must be consecutive"
118 call messages_fatal(1)
119 end if
120
121 select case (rdim)
122 case default
123 assert(.false.)
124 case (1)
125 call messages_not_implemented("The linear tetrahedron for 1D systems")
126 case (2)
127 block
128 integer, parameter :: submesh_triangles(2,3) = reshape([ &
129 1, 2, 3, &
130 1, 4, 3], shape(submesh_triangles), order=[2, 1])
131 ! coordinates of corners and neighbors in barycentric coordinates
132 integer, parameter :: b(10,3) = reshape([ &
133 1 , 0 , 0 , & ! k1
134 0 , 1 , 0 , & ! k2
135 0 , 0 , 1 , & ! k3
136 2 , -1, 0 , & ! 2 * k1 - k2
137 0 , 2 , -1, & ! 2 * k2 - k3
138 2 , 0 , -1, & ! 2 * k1 - k3
139 -1, 0 , 2 , & ! 2 * k3 - k1
140 -1, 2 , 0 , & ! 2 * k2 - k1
141 0 , -1, 2 , & ! 2 * k3 - k2
142 1 , -1, 1 & ! k1 - k2 + k3
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
148
149 this%n_simplices = 2 * npoints
150 this%dim = merge(10, 3, opt)
151 safe_allocate(this%simplices(this%n_simplices, this%dim))
152
153 do i = 1, raxis(1)
154 do j = 1, raxis(2)
155 ip1 = modulo(i, raxis(1)) + 1
156 jp1 = modulo(j, raxis(2)) + 1
157 corners(:,:) = reshape([ &
158 i , j , &
159 ip1 , j , &
160 ip1 , jp1 , &
161 i , jp1 ], shape(corners), order=[2, 1])
162
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), :)
169 do ik = 1, this%dim
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
174 end do
175 end do
176 end do
177 end do
178 end block
179 case(3)
180 block
181 integer, parameter :: submesh_tetras(6,4) = reshape([ &
182 1, 2, 3, 6, &
183 1, 3, 5, 6, &
184 3, 5, 6, 7, &
185 3, 6, 7, 8, &
186 3, 4, 6, 8, &
187 2, 3, 4, 6], shape(submesh_tetras), order=[2, 1])
188 ! coordinates of corners and neighbors in barycentric coordinates
189 integer, parameter :: b(20,4) = reshape([ &
190 1 , 0 , 0 , 0 , & ! k1
191 0 , 1 , 0 , 0 , & ! k2
192 0 , 0 , 1 , 0 , & ! k3
193 0 , 0 , 0 , 1 , & ! k4
194 2 , -1, 0 , 0 , & ! 2 * k1 - k2
195 0 , 2 , -1, 0 , & ! 2 * k2 - k3
196 0 , 0 , 2 , -1, & ! 2 * k3 - k4
197 -1, 0 , 0 , 2 , & ! 2 * k4 - k1
198 2 , 0 , -1, 0 , & ! 2 * k1 - k3
199 0 , 2 , 0 , -1, & ! 2 * k2 - k4
200 -1, 0 , 2 , 0 , & ! 2 * k3 - k1
201 0 , -1, 0 , 2 , & ! 2 * k4 - k2
202 2 , 0 , 0 , -1, & ! 2 * k1 - k4
203 -1, 2 , 0 , 0 , & ! 2 * k2 - k1
204 0 , -1, 2 , 0 , & ! 2 * k3 - k2
205 0 , 0 , -1, 2 , & ! 2 * k4 - k3
206 -1, 1 , 0 , 1 , & ! k4 - k1 + k2
207 1 , -1, 1 , 0 , & ! k1 - k2 + k3
208 0 , 1 , -1, 1 , & ! k2 - k3 + k4
209 1 , 0 , 1 , -1 & ! k3 - k4 + k1
210 ], shape(b), order=[2, 1])
211
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
215
216 this%n_simplices = 6 * npoints
217 this%dim = merge(20, 4, opt)
218 safe_allocate(this%simplices(this%n_simplices, this%dim))
219
220 do i = 1, raxis(1)
221 do j = 1, raxis(2)
222 do k = 1, raxis(3)
223 ip1 = modulo(i, raxis(1)) + 1
224 jp1 = modulo(j, raxis(2)) + 1
225 kp1 = modulo(k, raxis(3)) + 1
226 corners(:,:) = reshape([ &
227 i , j , k , &
228 ip1 , j , k , &
229 i , jp1 , k , &
230 ip1 , jp1 , k , &
231 i , j , kp1 , &
232 ip1 , j , kp1 , &
233 i , jp1 , kp1 , &
234 ip1 , jp1 , kp1 ], shape(corners), order=[2, 1])
235
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), :)
243 do ik = 1, this%dim
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
248 end do
249 end do
250 end do
251 end do
252 end do
253 end block
254 end select
255
256 safe_deallocate_a(kl123)
257
258 pop_sub(simplex_init)
259 end function simplex_init
260
264 subroutine simplex_end(this)
265 type(simplex_t), intent(inout) :: this
266 push_sub(simplex_end)
267 safe_deallocate_a(this%simplices)
268 pop_sub(simplex_end)
269 end subroutine simplex_end
270
280 subroutine simplex_weights_2d(etriangle, eF, weights, dos)
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
285
286 real(real64) :: e(3), e1, e2, e3
287 real(real64) :: w(3)
288 integer :: idx(3)
289
290 real(real64), parameter :: vt_vg = 1.0_real64 / 2.0_real64
291 real(real64), parameter :: vt_3vg = vt_vg / 3.0_real64
292
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])
298
299 ! no PUSH_SUB, called too often
300
301 select case (size(etriangle))
302 case (3)
303 e(:) = etriangle(:)
304 case (10)
305 e(:) = m_zero
306 block
307 integer :: i
308 do i = 1, size(etriangle)
309 e(:) = e(:) + p(:,i) * etriangle(i)
310 end do
311 end block
312 case default
313 assert(.false.)
314 end select
315
316 idx = [1,2,3]
317 call sort(e, idx)
318 e1 = e(1)
319 e2 = e(2)
320 e3 = e(3)
321
322 if (ef <= e1) then
323 w(:) = m_zero
324 dos = m_zero
325 elseif (e3 < ef) then
326 w(:) = vt_3vg
327 dos = m_zero
328 elseif (e1 < ef .and. ef <= e2) then
329 block
330 real(real64) :: e21, e31, c
331 e21 = e2 - e1
332 e31 = e3 - e1
333 c = vt_3vg * (ef - e1) ** 2 / (e21 * e31)
334
335 w(:) = c * [ &
336 3.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31), &
337 (ef - e1) / e21, &
338 (ef - e1) / e31]
339
340 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
341 end block
342 elseif (e2 < ef .and. ef <= e3) then
343 block
344 real(real64) :: e23, e31, c1, c2
345 e23 = e2 - e3
346 e31 = e3 - e1
347 c1 = vt_3vg
348 c2 = vt_3vg * (ef - e3) ** 2 / (e23 * e31)
349
350 w(:) = [ &
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))]
354
355 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
356 end block
357 else
358 assert(.false.)
359 end if
360
361 select case (size(etriangle))
362 case (3)
363 weights(idx) = w + dos * (sum(e) - 3.0_real64 * e) / 24.0_real64
364 case (10)
365 weights(idx) = w ! no Bloechl correction
366 case default
367 assert(.false.)
368 end select
369 end subroutine simplex_weights_2d
370
376 subroutine simplex_dos_2d(etriangle, eF, dos)
377 real(real64), intent(in) :: etriangle(:)
378 real(real64), intent(in) :: ef
379 real(real64), intent(out) :: dos
380
381 real(real64) :: e(3), e1, e2, e3
382
383 real(real64), parameter :: vt_vg = 1.0_real64 / 2.0_real64
384
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])
390
391 ! no PUSH_SUB, called too often
392
393 select case (size(etriangle))
394 case (3)
395 e(:) = etriangle(:)
396 case (10)
397 e(:) = m_zero
398 block
399 integer :: i
400 do i = 1, size(etriangle)
401 e(:) = e(:) + p(:,i) * etriangle(i)
402 end do
403 end block
404 case default
405 assert(.false.)
406 end select
407
408 call sort(e)
409 e1 = e(1)
410 e2 = e(2)
411 e3 = e(3)
412
413 if (ef <= e1 .or. e3 < ef) then
414 dos = m_zero
415 elseif (e1 < ef .and. ef <= e2) then
416 block
417 real(real64) :: e21, e31
418 e21 = e2 - e1
419 e31 = e3 - e1
420 dos = 2.0_real64 * vt_vg * (ef - e1) / (e21 * e31)
421 end block
422 elseif (e2 < ef .and. ef <= e3) then
423 block
424 real(real64) :: e23, e31
425 e23 = e2 - e3
426 e31 = e3 - e1
427 dos = 2.0_real64 * vt_vg * (ef - e3) / (e23 * e31)
428 end block
429 else
430 assert(.false.)
431 end if
432 end subroutine simplex_dos_2d
433
446 subroutine simplex_weights_3d(etetra, eF, weights, dos)
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
451
452 real(real64) :: e(4), e1, e2, e3, e4
453 real(real64) :: w(4)
454 integer :: idx(4)
455
456 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
457 real(real64), parameter :: vt_4vg = vt_vg / 4.0_real64
458
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])
465
466 ! no PUSH_SUB, called too often
467
468 select case (size(etetra))
469 case (4)
470 e(:) = etetra(:)
471 case (20)
472 e(:) = m_zero
473 block
474 integer :: i
475 do i = 1, size(etetra)
476 e(:) = e(:) + p(:,i) * etetra(i)
477 end do
478 end block
479 case default
480 assert(.false.)
481 end select
482
483 idx = [1,2,3,4]
484 call sort(e, idx)
485 e1 = e(1)
486 e2 = e(2)
487 e3 = e(3)
488 e4 = e(4)
489
490 if (e1 >= ef) then
491 w(:) = m_zero
492 dos = m_zero
493 elseif (e4 < ef) then
494 w(:) = vt_4vg
495 dos = m_zero
496 elseif (e1 < ef .and. ef <= e2) then
497 block
498 real(real64) :: e21, e31, e41, c
499 e21 = e2 - e1
500 e31 = e3 - e1
501 e41 = e4 - e1
502 c = vt_4vg * (ef - e1) ** 3 / (e21 * e31 * e41)
503
504 w(:) = c * [ &
505 4.0_real64 - (ef - e1) * (m_one / e21 + m_one / e31 + m_one / e41), &
506 (ef - e1) / e21, &
507 (ef - e1) / e31, &
508 (ef - e1) / e41]
509
510 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
511 end block
512 elseif (e2 < ef .and. ef <= e3) then
513 block
514 real(real64) :: e21, e31, e32, e41, e42, c1, c2, c3
515 e21 = e2 - e1
516 e31 = e3 - e1
517 e32 = e3 - e2
518 e41 = e4 - e1
519 e42 = e4 - e2
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)
523
524 w(:) = [ &
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]
529
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))
532 end block
533 elseif (e3 < ef .and. ef <= e4) then
534 block
535 real(real64) :: e41, e42, e43, c
536 e41 = e4 - e1
537 e42 = e4 - e2
538 e43 = e4 - e3
539 c = vt_4vg * (e4 - ef) ** 3 / (e41 * e42 * e43)
540
541 w(:) = vt_4vg - c * [ &
542 (e4 - ef) / e41, &
543 (e4 - ef) / e42, &
544 (e4 - ef) / e43, &
545 4.0_real64 - (e4 - ef) * (m_one / e41 + m_one / e42 + m_one / e43)]
546
547 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
548 end block
549 else
550 assert(.false.)
551 end if
552
553 select case (size(etetra))
554 case (4)
555 weights(idx) = w + dos * (sum(e) - 4.0_real64 * e) / 40.0_real64
556 case (20)
557 weights(idx) = w ! no Bloechl correction
558 case default
559 assert(.false.)
560 end select
561 end subroutine simplex_weights_3d
562
574 subroutine simplex_dos_3d(etetra, eF, dos)
575 real(real64), intent(in) :: etetra(:)
576 real(real64), intent(in) :: ef
577 real(real64), intent(out) :: dos
578
579 real(real64) :: e(4), e1, e2, e3, e4
580
581 real(real64), parameter :: vt_vg = 1.0_real64 / 6.0_real64
582
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])
589
590 ! no PUSH_SUB, called too often
591
592 select case (size(etetra))
593 case (4)
594 e(:) = etetra(:)
595 case (20)
596 e(:) = m_zero
597 block
598 integer :: i
599 do i = 1, size(etetra)
600 e(:) = e(:) + p(:,i) * etetra(i)
601 end do
602 end block
603 case default
604 assert(.false.)
605 end select
606
607 call sort(e)
608 e1 = e(1)
609 e2 = e(2)
610 e3 = e(3)
611 e4 = e(4)
612
613 if (e1 >= ef .or. e4 < ef) then
614 dos = m_zero
615 elseif (e1 < ef .and. ef <= e2) then
616 block
617 real(real64) :: e21, e31, e41
618 e21 = e2 - e1
619 e31 = e3 - e1
620 e41 = e4 - e1
621 dos = vt_vg * 3.0_real64 * (ef - e1) ** 2 / (e21 * e31 * e41)
622 end block
623 elseif (e2 < ef .and. ef <= e3) then
624 block
625 real(real64) :: e21, e31, e32, e41, e42
626 e21 = e2 - e1
627 e31 = e3 - e1
628 e32 = e3 - e2
629 e41 = e4 - e1
630 e42 = e4 - e2
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))
633 end block
634 elseif (e3 < ef .and. ef <= e4) then
635 block
636 real(real64) :: e41, e42, e43
637 e41 = e4 - e1
638 e42 = e4 - e2
639 e43 = e4 - e3
640 dos = vt_vg * 3.0_real64 * (e4 - ef) ** 2 / (e41 * e42 * e43)
641 end block
642 else
643 assert(.false.)
644 end if
645 end subroutine simplex_dos_3d
646
647end 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_2d(etriangle, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
Definition: simplex.F90:376
subroutine, public simplex_dos_3d(etetra, eF, dos)
Get only the DOS contribution of a single tetrahedron.
Definition: simplex.F90:670
type(simplex_t) function, pointer, public simplex_init(dim, naxis, nshifts, shift, kpoints, equiv, opt)
Constructor for linear simplex methods.
Definition: simplex.F90:172
subroutine, public simplex_end(this)
Destructor for linear simplex methods.
Definition: simplex.F90:360
subroutine, public simplex_dos_2d(etriangle, eF, dos)
Get only the DOS contribution of a single triangle.
Definition: simplex.F90:472
subroutine, public simplex_weights_3d(etetra, eF, weights, dos)
Get the weights and DOS contribution of a single tetrahedron.
Definition: simplex.F90:542
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119