The Elk Code
Loading...
Searching...
No Matches
genffacgp.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2008 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: genffacgp
8
! !INTERFACE:
9
pure
subroutine
genffacgp
(ngp,gpc,ld,ffacgp)
10
! !USES:
11
use
modmain
12
! !INPUT/OUTPUT PARAMETERS:
13
! ngp : number of G+p-vectors (in,integer)
14
! gpc : length of G+p-vectors (in,real(ngp))
15
! ld : leading dimension (in,integer)
16
! ffacgp : form factors (out,real(ld,nspecies))
17
! !DESCRIPTION:
18
! Generates the form factors used to determine the smooth characteristic
19
! function. See {\tt gencfun} for details.
20
!
21
! !REVISION HISTORY:
22
! Created January 2003 (JKD)
23
!EOP
24
!BOC
25
implicit none
26
! arguments
27
integer
,
intent(in)
:: ngp
28
real
(8),
intent(in)
:: gpc(ngp)
29
integer
,
intent(in)
:: ld
30
real
(8),
intent(out)
:: ffacgp(ld,
nspecies
)
31
! local variables
32
integer
is,ig
33
real
(8) r,t0,t1
34
t0=
fourpi
/
omega
35
do
is=1,
nspecies
36
r=
rmt
(is)
37
do
ig=1,ngp
38
if
(gpc(ig) >
epslat
)
then
39
t1=gpc(ig)*r
40
ffacgp(ig,is)=t0*(sin(t1)-t1*cos(t1))/(gpc(ig)**3)
41
else
42
ffacgp(ig,is)=(t0/3.d0)*(r**3)
43
end if
44
end do
45
end do
46
end subroutine
47
!EOC
48
genffacgp
pure subroutine genffacgp(ngp, gpc, ld, ffacgp)
Definition
genffacgp.f90:10
modmain
Definition
modmain.f90:6
modmain::rmt
real(8), dimension(maxspecies) rmt
Definition
modmain.f90:162
modmain::omega
real(8) omega
Definition
modmain.f90:20
modmain::fourpi
real(8), parameter fourpi
Definition
modmain.f90:1231
modmain::epslat
real(8) epslat
Definition
modmain.f90:24
modmain::nspecies
integer nspecies
Definition
modmain.f90:34
genffacgp.f90
Generated by
1.9.8