The Elk Code
rhomagv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
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 rhomagv
7 use modmain
8 use modmpi
9 use modomp
10 implicit none
11 ! local variables
12 integer ik,ispn,idm
13 integer is,ias,n,nthd
14 ! allocatable arrays
15 real(8), allocatable :: rhomt_(:,:),rhoir_(:),magmt_(:,:,:),magir_(:,:)
16 complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
17 ! set the charge density and magnetisation to zero
18 allocate(rhomt_(npcmtmax,natmtot),rhoir_(ngtc))
19 rhomt_(:,:)=0.d0
20 rhoir_(:)=0.d0
21 if (spinpol) then
22  allocate(magmt_(npcmtmax,natmtot,ndmag),magir_(ngtc,ndmag))
23 else
24  allocate(magmt_(1,1,1),magir_(1,1))
25 end if
26 magmt_(:,:,:)=0.d0
27 magir_(:,:)=0.d0
28 call holdthd(nkpt/np_mpi,nthd)
29 !$OMP PARALLEL DEFAULT(SHARED) &
30 !$OMP PRIVATE(apwalm,evecfv,evecsv,ispn) &
31 !$OMP REDUCTION(+:rhomt_,rhoir_,magmt_,magir_) &
32 !$OMP NUM_THREADS(nthd)
33 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
34 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
35 !$OMP DO SCHEDULE(STATIC)
36 do ik=1,nkpt
37 ! distribute among MPI processes
38  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
39 ! get the eigenvectors from file
40  call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
41  call getevecsv(filext,ik,vkl(:,ik),evecsv)
42 ! find the matching coefficients
43  do ispn=1,nspnfv
44  call match(ngk(ispn,ik),vgkc(:,:,ispn,ik),gkc(:,ispn,ik), &
45  sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
46  end do
47 ! add to the density and magnetisation
48  call rhomagk(ngk(:,ik),igkig(:,:,ik),wkpt(ik),occsv(:,ik),apwalm,evecfv, &
49  evecsv,rhomt_,rhoir_,magmt_,magir_)
50 end do
51 !$OMP END DO
52 deallocate(apwalm,evecfv,evecsv)
53 !$OMP END PARALLEL
54 call freethd(nthd)
55 ! add density from each process and redistribute
56 if (np_mpi > 1) then
58  call mpi_allreduce(mpi_in_place,rhomt_,n,mpi_double_precision,mpi_sum,mpicom,&
59  ierror)
60  call mpi_allreduce(mpi_in_place,rhoir_,ngtc,mpi_double_precision,mpi_sum, &
61  mpicom,ierror)
62  if (spinpol) then
63  n=n*ndmag
64  call mpi_allreduce(mpi_in_place,magmt_,n,mpi_double_precision,mpi_sum, &
65  mpicom,ierror)
66  n=ngtc*ndmag
67  call mpi_allreduce(mpi_in_place,magir_,n,mpi_double_precision,mpi_sum, &
68  mpicom,ierror)
69  end if
70 end if
71 ! copy density to the global arrays
72 do ias=1,natmtot
73  is=idxis(ias)
74  rhomt(1:npcmt(is),ias)=rhomt_(1:npcmt(is),ias)
75 end do
76 rhoir(1:ngtc)=rhoir_(1:ngtc)
77 if (spinpol) then
78  do idm=1,ndmag
79  do ias=1,natmtot
80  is=idxis(ias)
81  magmt(1:npcmt(is),ias,idm)=magmt_(1:npcmt(is),ias,idm)
82  end do
83  magir(1:ngtc,idm)=magir_(1:ngtc,idm)
84  end do
85 end if
86 deallocate(rhomt_,rhoir_,magmt_,magir_)
87 ! convert muffin-tin density/magnetisation to spherical harmonics
88 call rhomagsh
89 call holdthd(2,nthd)
90 !$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
91 !$OMP PRIVATE(idm) &
92 !$OMP NUM_THREADS(nthd)
93 !$OMP SECTION
94 ! symmetrise the density
96  rhoir)
97 ! convert the muffin-tin density from coarse to fine radial mesh
98 call rfmtctof(rhomt)
99 ! convert the interstitial density from coarse to fine grid
100 call rfirctof(rhoir,rhoir)
101 !$OMP SECTION
102 if (spinpol) then
103 ! symmetrise the magnetisation
106 ! convert the muffin-tin magnetisation from coarse to fine radial mesh
107  do idm=1,ndmag
108  call rfmtctof(magmt(:,:,idm))
109  end do
110 ! convert the interstitial magnetisation from coarse to fine grid
111  do idm=1,ndmag
112  call rfirctof(magir(:,idm),magir(:,idm))
113  end do
114 end if
115 !$OMP END PARALLEL SECTIONS
116 call freethd(nthd)
117 end subroutine
118 
integer nmatmax
Definition: modmain.f90:853
subroutine rhomagsh
Definition: rhomagsh.f90:10
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
subroutine rfirctof(rfirc, rfir)
Definition: rfirctof.f90:10
integer npcmtmax
Definition: modmain.f90:216
character(256) filext
Definition: modmain.f90:1301
integer nfgrzc
Definition: modmain.f90:414
integer ngtot
Definition: modmain.f90:390
integer ngtc
Definition: modmain.f90:392
logical spinpol
Definition: modmain.f90:228
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
integer lmmaxapw
Definition: modmain.f90:199
integer nkpt
Definition: modmain.f90:461
integer ndmag
Definition: modmain.f90:238
integer ngkmax
Definition: modmain.f90:499
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
Definition: modomp.f90:6
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
subroutine symrf(nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld, rfmt, rfir)
Definition: symrf.f90:11
integer, dimension(:), allocatable igrzfc
Definition: modmain.f90:418
subroutine rhomagv
Definition: rhomagv.f90:7
integer ngvc
Definition: modmain.f90:398
subroutine rhomagk(ngp, igpig, wppt, occsvp, apwalm, evecfv, evecsv, rhomt_, rhoir_, magmt_, magir_)
Definition: rhomagk.f90:11
subroutine rfmtctof(rfmt)
Definition: rfmtctof.f90:10
integer np_mpi
Definition: modmpi.f90:13
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
subroutine symrvf(tspin, tnc, nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld1, rvfmt, ld2, rvfir)
Definition: symrvf.f90:11
integer nstsv
Definition: modmain.f90:889
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
real(8), dimension(:), allocatable wkpt
Definition: modmain.f90:475
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer apwordmax
Definition: modmain.f90:760
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
Definition: modmpi.f90:6
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer, dimension(3) ngdgc
Definition: modmain.f90:388
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer npmtmax
Definition: modmain.f90:216
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
logical ncmag
Definition: modmain.f90:240
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nstfv
Definition: modmain.f90:887
integer mpicom
Definition: modmpi.f90:11
integer nspnfv
Definition: modmain.f90:289
integer ierror
Definition: modmpi.f90:19