The Elk Code
Loading...
Searching...
No Matches
gengclgq.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2017 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
pure
subroutine
gengclgq
(treg,iq,ngq,gqc,gclgq)
7
use
modmain
8
implicit none
9
! arguments
10
logical
,
intent(in)
:: treg
11
integer
,
intent(in)
:: iq,ngq
12
real
(8),
intent(in)
:: gqc(ngq)
13
real
(8),
intent(out)
:: gclgq(ngq)
14
! local variables
15
integer
ig
16
real
(8) t1,t2
17
if
(treg)
then
18
! regularise 1/(G+q)² for G+q in the first Brillouin zone
19
t1=sqrt(
vqc
(1,iq)**2+
vqc
(2,iq)**2+
vqc
(3,iq)**2)
20
do
ig=1,ngq
21
t2=gqc(ig)
22
if
(abs(t1-t2) <
epslat
)
then
23
gclgq(ig)=
gclq
(iq)
24
else
if
(t2 >
epslat
)
then
25
gclgq(ig)=
fourpi
/t2**2
26
else
27
gclgq(ig)=0.d0
28
end if
29
end do
30
else
31
! no regularisation
32
do
ig=1,ngq
33
t1=gqc(ig)
34
if
(t1 >
epslat
)
then
35
gclgq(ig)=
fourpi
/t1**2
36
else
37
gclgq(ig)=0.d0
38
end if
39
end do
40
end if
41
end subroutine
42
gengclgq
pure subroutine gengclgq(treg, iq, ngq, gqc, gclgq)
Definition
gengclgq.f90:7
modmain
Definition
modmain.f90:6
modmain::gclq
real(8), dimension(:), allocatable gclq
Definition
modmain.f90:553
modmain::fourpi
real(8), parameter fourpi
Definition
modmain.f90:1231
modmain::vqc
real(8), dimension(:,:), allocatable vqc
Definition
modmain.f90:547
modmain::epslat
real(8) epslat
Definition
modmain.f90:24
gengclgq.f90
Generated by
1.9.8