The Elk Code
Loading...
Searching...
No Matches
writegclq.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 General Public License.
4
! See the file COPYING for license details.
5
6
!BOP
7
! !ROUTINE: writegclq
8
! !INTERFACE:
9
subroutine
writegclq
10
! !USES:
11
use
modmain
12
! !DESCRIPTION:
13
! Outputs the volume-averaged integral of $4\pi/q^2$ in the small
14
! parallelepiped around each discrete $q$-point to the file {\tt GCLQ.OUT}.
15
! These represent the regularised Coulomb Green's function in reciprocal
16
! space for small $q$. See the routine gengclq.
17
!
18
! !REVISION HISTORY:
19
! Created June 2005 (JKD)
20
!EOP
21
!BOC
22
implicit none
23
! local variables
24
integer
iq
25
real
(8) t1
26
open
(50,file=
'GCLQ'
//trim(
filext
),form=
'FORMATTED'
,action=
'WRITE'
)
27
write
(50,π²
'(I6," : nqpt; q-point, vql, gclq, 4/q below")'
)
nqpt
28
do
iq=1,
nqpt
29
t1=
vqc
(1,iq)**2+
vqc
(2,iq)**2+
vqc
(3,iq)**2
30
if
(t1 > 1.d-12)
then
31
t1=
fourpi
/t1
32
else
33
t1=0.d0
34
end if
35
write
(50,
'(I6,5G18.10)'
) iq,
vql
(:,iq),
gclq
(iq),t1
36
end do
37
close
(50)
38
end subroutine
39
!EOC
40
modmain
Definition
modmain.f90:6
modmain::gclq
real(8), dimension(:), allocatable gclq
Definition
modmain.f90:553
modmain::filext
character(256) filext
Definition
modmain.f90:1300
modmain::fourpi
real(8), parameter fourpi
Definition
modmain.f90:1231
modmain::vqc
real(8), dimension(:,:), allocatable vqc
Definition
modmain.f90:547
modmain::nqpt
integer nqpt
Definition
modmain.f90:525
modmain::vql
real(8), dimension(:,:), allocatable vql
Definition
modmain.f90:545
writegclq
subroutine writegclq
Definition
writegclq.f90:10
writegclq.f90
Generated by
1.9.8