The Elk Code
rhocoreu.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 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 
6 subroutine rhocoreu
7 use modmain
8 use modulr
9 use modomp
10 implicit none
11 ! local variables
12 integer is,ias,ir,i
13 integer nrc,nrci,irc
14 integer npc,n,nthd
15 real(8) t1
16 ! allocatable arrays
17 real(8), allocatable :: rfmt(:,:)
18 ! generate the core density in spherical coordinates
19 allocate(rfmt(npcmtmax,natmtot))
20 do ias=1,natmtot
21  is=idxis(ias)
22  nrc=nrcmt(is)
23  nrci=nrcmti(is)
24  n=lmmaxi-1
25  ir=1
26  i=1
27  do irc=1,nrci
28  t1=rhocr(ir,ias,1)*y00
29  rfmt(i:i+n,ias)=t1
30  ir=ir+lradstp
31  i=i+lmmaxi
32  end do
33  n=lmmaxo-1
34  do irc=nrci+1,nrc
35  t1=rhocr(ir,ias,1)*y00
36  rfmt(i:i+n,ias)=t1
37  ir=ir+lradstp
38  i=i+lmmaxo
39  end do
40 end do
41 ! add to the ultra long-range density
42 call holdthd(nqpt,nthd)
43 !$OMP PARALLEL DO DEFAULT(SHARED) &
44 !$OMP PRIVATE(ias,is,npc) &
45 !$OMP SCHEDULE(DYNAMIC) &
46 !$OMP NUM_THREADS(nthd)
47 do ir=1,nqpt
48  do ias=1,natmtot
49  is=idxis(ias)
50  npc=npcmt(is)
51  rhormt(1:npc,ias,ir)=rhormt(1:npc,ias,ir)+rfmt(1:npc,ias)
52  end do
53 end do
54 !$OMP END PARALLEL DO
55 call freethd(nthd)
56 deallocate(rfmt)
57 end subroutine
58 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer npcmtmax
Definition: modmain.f90:216
integer lmmaxo
Definition: modmain.f90:203
integer nqpt
Definition: modmain.f90:525
Definition: modomp.f90:6
real(8), dimension(:,:,:), allocatable rhocr
Definition: modmain.f90:941
subroutine rhocoreu
Definition: rhocoreu.f90:7
integer lradstp
Definition: modmain.f90:171
real(8), dimension(:,:,:), pointer, contiguous rhormt
Definition: modulr.f90:53
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
Definition: modulr.f90:6
real(8), parameter y00
Definition: modmain.f90:1236