The Elk Code
 
Loading...
Searching...
No Matches
genjlgpr.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
6subroutine genjlgpr(ngp,gpc,jlgpr)
7use modmain
8use modomp
9implicit none
10! arguments
11integer, intent(in) :: ngp
12real(8), intent(in) :: gpc(ngp)
13real(8), intent(out) :: jlgpr(njcmax,nspecies,ngp)
14! local variables
15integer ig,is,n,i,nthd
16integer nrc,nrci,irc
17real(8) t1,t2
18! generate spherical Bessel functions on the coarse radial mesh over all species
19call holdthd(ngp,nthd)
20!$OMP PARALLEL DO DEFAULT(SHARED) &
21!$OMP PRIVATE(t1,t2,is,nrc) &
22!$OMP PRIVATE(nrci,n,i,irc) &
23!$OMP SCHEDULE(DYNAMIC) &
24!$OMP NUM_THREADS(nthd)
25do ig=1,ngp
26 t1=gpc(ig)
27 do is=1,nspecies
28 nrc=nrcmt(is)
29 nrci=nrcmti(is)
30 n=lmaxi+1
31 i=1
32 do irc=1,nrci
33 t2=t1*rcmt(irc,is)
34 call sbessel(lmaxi,t2,jlgpr(i,is,ig))
35 i=i+n
36 end do
37 n=lmaxo+1
38 do irc=nrci+1,nrc
39 t2=t1*rcmt(irc,is)
40 call sbessel(lmaxo,t2,jlgpr(i,is,ig))
41 i=i+n
42 end do
43 end do
44end do
45!$OMP END PARALLEL DO
46call freethd(nthd)
47end subroutine
48
subroutine genjlgpr(ngp, gpc, jlgpr)
Definition genjlgpr.f90:7
real(8), dimension(:,:), allocatable rcmt
Definition modmain.f90:177
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer lmaxo
Definition modmain.f90:201
integer lmaxi
Definition modmain.f90:205
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine sbessel(lmax, x, jl)
Definition sbessel.f90:10