The Elk Code
genjlgprmt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2006 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: genjlgprmt
8 ! !INTERFACE:
9 subroutine genjlgprmt(lmax,ngp,gpc,ld,jlgprmt)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! lmax : angular momentum cut-off (in,integer)
14 ! ngp : number of G+p-vectors (in,integer)
15 ! gpc : length of G+p-vectors (in,real(ngkmax))
16 ! ld : leading dimension (in,integer)
17 ! jlgprmt : spherical Bessel functions (out,real(0:lmax,ld,nspecies))
18 ! !DESCRIPTION:
19 ! Calculates and stores the spherical Bessel functions
20 ! $j_l(|{\bf G}+{\bf p}|{\bf R}_{\rm MT})$ for all input ${\bf G}+{\bf p}$
21 ! vectors and the muffin-tin radii ${\bf R}_{\rm MT}$ of every atomic species.
22 !
23 ! !REVISION HISTORY:
24 ! Created April 2002 (JKD)
25 !EOP
26 !BOC
27 implicit none
28 ! arguments
29 integer, intent(in) :: lmax,ngp
30 real(8), intent(in) :: gpc(ngp)
31 integer, intent(in) :: ld
32 real(8), intent(out) :: jlgprmt(0:lmax,ld,nspecies)
33 ! local variables
34 integer is,ig
35 real(8) r,x
36 do is=1,nspecies
37  r=rmt(is)
38  do ig=1,ngp
39  x=gpc(ig)*r
40  call sbessel(lmax,x,jlgprmt(:,ig,is))
41  end do
42 end do
43 end subroutine
44 !EOC
45 
subroutine sbessel(lmax, x, jl)
Definition: sbessel.f90:10
real(8), dimension(maxspecies) rmt
Definition: modmain.f90:162
subroutine genjlgprmt(lmax, ngp, gpc, ld, jlgprmt)
Definition: genjlgprmt.f90:10