The Elk Code
 
Loading...
Searching...
No Matches
brzint.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: brzint
8! !INTERFACE:
9subroutine brzint(nsm,ngridk,nsk,ivkik,nw,wint,n,ld,e,f,g)
10! !USES:
11use modomp
12! !INPUT/OUTPUT PARAMETERS:
13! nsm : level of smoothing for output function (in,integer)
14! ngridk : k-point grid size (in,integer(3))
15! nsk : k-point subdivision grid size (in,integer(3))
16! ivkik : map from (i1,i2,i3) to k-point index
17! (in,integer(0:ngridk(1)-1,0:ngridk(2)-1,0:ngridk(3)-1))
18! nw : number of energy divisions (in,integer)
19! wint : energy interval (in,real(2))
20! n : number of functions to integrate (in,integer)
21! ld : leading dimension (in,integer)
22! e : array of energies as a function of k-points (in,real(ld,*))
23! f : array of weights as a function of k-points (in,real(ld,*))
24! g : output function (out,real(nw))
25! !DESCRIPTION:
26! Given energy and weight functions, $e$ and $f$, on the Brillouin zone and a
27! set of equidistant energies $\omega_i$, this routine computes the integrals
28! $$ g(\omega_i)=\frac{\Omega}{(2\pi)^3}\int_{\rm BZ} f({\bf k})
29! \delta(\omega_i-e({\bf k}))d{\bf k}, $$
30! where $\Omega$ is the unit cell volume. This is done by first interpolating
31! $e$ and $f$ on a finer $k$-point grid using the trilinear method. Then for
32! each $e({\bf k})$ on the finer grid the nearest $\omega_i$ is found and
33! $f({\bf k})$ is accumulated in $g(\omega_i)$. If the output function is
34! noisy then either {\tt nsk} should be increased or {\tt nw} decreased.
35! Alternatively, the output function can be artificially smoothed up to a
36! level given by {\tt nsm}. See routine {\tt fsmooth}.
37!
38! !REVISION HISTORY:
39! Created October 2003 (JKD)
40! Improved efficiency, May 2007 (Sebastian Lebegue)
41! Added parallelism, March 2020 (JKD)
42!EOP
43!BOC
44implicit none
45! arguments
46integer, intent(in) :: nsm,ngridk(3),nsk(3)
47integer, intent(in) :: ivkik(0:ngridk(1)-1,0:ngridk(2)-1,0:ngridk(3)-1)
48integer, intent(in) :: nw
49real(8), intent(in) :: wint(2)
50integer, intent(in) :: n,ld
51real(8), intent(in) :: e(ld,*),f(ld,*)
52real(8), intent(out) :: g(nw)
53! local variables
54integer nk,i1,i2,i3,j1,j2,j3,k1,k2,k3,i,iw,nthd
55integer i000,i001,i010,i011,i100,i101,i110,i111
56real(8) wd,dw,dwi,w1,t1,t2
57! automatic arrays
58real(8) f0(n),f1(n),e0(n),e1(n)
59real(8) f00(n),f01(n),f10(n),f11(n)
60real(8) e00(n),e01(n),e10(n),e11(n)
61if ((ngridk(1) < 1).or.(ngridk(2) < 1).or.(ngridk(3) < 1)) then
62 write(*,*)
63 write(*,'("Error(brzint): ngridk < 1 : ",3I8)') ngridk
64 write(*,*)
65 stop
66end if
67if ((nsk(1) < 1).or.(nsk(2) < 1).or.(nsk(3) < 1)) then
68 write(*,*)
69 write(*,'("Error(brzint): nsk < 1 : ",3I8)') nsk
70 write(*,*)
71 stop
72end if
73! total number of k-points
74nk=ngridk(1)*ngridk(2)*ngridk(3)
75! length of interval
76wd=wint(2)-wint(1)
77! energy step size
78dw=wd/dble(nw)
79dwi=1.d0/dw
80w1=wint(1)*dwi
81g(1:nw)=0.d0
82call holdthd(nk,nthd)
83!$OMP PARALLEL DO DEFAULT(SHARED) &
84!$OMP PRIVATE(f0,f1,e0,e1) &
85!$OMP PRIVATE(f00,f01,f10,f11) &
86!$OMP PRIVATE(e00,e01,e10,e11) &
87!$OMP PRIVATE(k1,k2,k3,i1,i2,i3) &
88!$OMP PRIVATE(i000,i001,i010,i011) &
89!$OMP PRIVATE(i100,i101,i110,i111) &
90!$OMP PRIVATE(t1,t2,i,iw) &
91!$OMP REDUCTION(+:g) &
92!$OMP SCHEDULE(DYNAMIC) &
93!$OMP NUM_THREADS(nthd) &
94!$OMP COLLAPSE(3)
95do j1=0,ngridk(1)-1
96 do j2=0,ngridk(2)-1
97 do j3=0,ngridk(3)-1
98 k1=mod(j1+1,ngridk(1))
99 k2=mod(j2+1,ngridk(2))
100 k3=mod(j3+1,ngridk(3))
101 i000=ivkik(j1,j2,j3); i001=ivkik(j1,j2,k3)
102 i010=ivkik(j1,k2,j3); i011=ivkik(j1,k2,k3)
103 i100=ivkik(k1,j2,j3); i101=ivkik(k1,j2,k3)
104 i110=ivkik(k1,k2,j3); i111=ivkik(k1,k2,k3)
105 do i1=0,nsk(1)-1
106 t2=dble(i1)/dble(nsk(1))
107 t1=1.d0-t2
108 f00(1:n)=f(1:n,i000)*t1+f(1:n,i100)*t2
109 f01(1:n)=f(1:n,i001)*t1+f(1:n,i101)*t2
110 f10(1:n)=f(1:n,i010)*t1+f(1:n,i110)*t2
111 f11(1:n)=f(1:n,i011)*t1+f(1:n,i111)*t2
112 t1=t1*dwi
113 t2=t2*dwi
114 e00(1:n)=e(1:n,i000)*t1+e(1:n,i100)*t2-w1
115 e01(1:n)=e(1:n,i001)*t1+e(1:n,i101)*t2-w1
116 e10(1:n)=e(1:n,i010)*t1+e(1:n,i110)*t2-w1
117 e11(1:n)=e(1:n,i011)*t1+e(1:n,i111)*t2-w1
118 do i2=0,nsk(2)-1
119 t2=dble(i2)/dble(nsk(2))
120 t1=1.d0-t2
121 f0(1:n)=f00(1:n)*t1+f10(1:n)*t2
122 f1(1:n)=f01(1:n)*t1+f11(1:n)*t2
123 e0(1:n)=e00(1:n)*t1+e10(1:n)*t2
124 e1(1:n)=e01(1:n)*t1+e11(1:n)*t2
125 do i3=0,nsk(3)-1
126 t2=dble(i3)/dble(nsk(3))
127 t1=1.d0-t2
128 do i=1,n
129 iw=nint(e0(i)*t1+e1(i)*t2)+1
130 if ((iw >= 1).and.(iw <= nw)) g(iw)=g(iw)+f0(i)*t1+f1(i)*t2
131 end do
132 end do
133 end do
134 end do
135 end do
136 end do
137end do
138!$OMP END PARALLEL DO
139call freethd(nthd)
140! normalise function
141t1=dw*dble(nk)*dble(nsk(1)*nsk(2)*nsk(3))
142t1=1.d0/t1
143g(1:nw)=t1*g(1:nw)
144! smooth output function if required
145if (nsm > 0) call fsmooth(nsm,nw,g)
146end subroutine
147!EOC
148
subroutine brzint(nsm, ngridk, nsk, ivkik, nw, wint, n, ld, e, f, g)
Definition brzint.f90:10
pure subroutine fsmooth(m, n, f)
Definition fsmooth.f90:10
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106