The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine rhocoreu
7use modmain
8use modulr
9use modomp
10implicit none
11! local variables
12integer is,ias,ir,i
13integer nrc,nrci,irc
14integer npc,n,nthd
15real(8) t1
16! allocatable arrays
17real(8), allocatable :: rfmt(:,:)
18! generate the core density in spherical coordinates
19allocate(rfmt(npcmtmax,natmtot))
20do 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
40end do
41! add to the ultra long-range density
42call holdthd(nqpt,nthd)
43!$OMP PARALLEL DO DEFAULT(SHARED) &
44!$OMP PRIVATE(ias,is,npc) &
45!$OMP SCHEDULE(DYNAMIC) &
46!$OMP NUM_THREADS(nthd)
47do 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
53end do
54!$OMP END PARALLEL DO
55call freethd(nthd)
56deallocate(rfmt)
57end subroutine
58
real(8), parameter y00
Definition modmain.f90:1233
integer lmmaxi
Definition modmain.f90:207
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer lradstp
Definition modmain.f90:171
integer nqpt
Definition modmain.f90:525
integer natmtot
Definition modmain.f90:40
integer npcmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:), allocatable rhocr
Definition modmain.f90:938
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
real(8), dimension(:,:,:), allocatable rhormt
Definition modulr.f90:52
subroutine rhocoreu
Definition rhocoreu.f90:7