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:
9pure subroutine genffacgp(ngp,gpc,ld,ffacgp)
10! !USES:
11use 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
25implicit none
26! arguments
27integer, intent(in) :: ngp
28real(8), intent(in) :: gpc(ngp)
29integer, intent(in) :: ld
30real(8), intent(out) :: ffacgp(ld,nspecies)
31! local variables
32integer is,ig
33real(8) r,t0,t1
35do 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
45end do
46end subroutine
47!EOC
48
pure subroutine genffacgp(ngp, gpc, ld, ffacgp)
Definition genffacgp.f90:10
real(8), dimension(maxspecies) rmt
Definition modmain.f90:162
real(8) omega
Definition modmain.f90:20
real(8), parameter fourpi
Definition modmain.f90:1231
real(8) epslat
Definition modmain.f90:24
integer nspecies
Definition modmain.f90:34