The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine rhomagv
7use modmain
8use modmpi
9use modomp
10implicit none
11! local variables
12integer ik,ispn,idm
13integer is,ias,n,nthd
14! allocatable arrays
15real(8), allocatable :: rhomt_(:,:),rhoir_(:),magmt_(:,:,:),magir_(:,:)
16complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
17! set the charge density and magnetisation to zero
18allocate(rhomt_(npcmtmax,natmtot),rhoir_(ngtc))
19rhomt_(:,:)=0.d0
20rhoir_(:)=0.d0
21if (spinpol) then
22 allocate(magmt_(npcmtmax,natmtot,ndmag),magir_(ngtc,ndmag))
23else
24 allocate(magmt_(1,1,1),magir_(1,1))
25end if
26magmt_(:,:,:)=0.d0
27magir_(:,:)=0.d0
28call 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)
33allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
34allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
35!$OMP DO SCHEDULE(STATIC)
36do 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_)
50end do
51!$OMP END DO
52deallocate(apwalm,evecfv,evecsv)
53!$OMP END PARALLEL
54call freethd(nthd)
55! add density from each process and redistribute
56if (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, &
62 if (spinpol) then
63 n=n*ndmag
64 call mpi_allreduce(mpi_in_place,magmt_,n,mpi_double_precision,mpi_sum, &
66 n=ngtc*ndmag
67 call mpi_allreduce(mpi_in_place,magir_,n,mpi_double_precision,mpi_sum, &
69 end if
70end if
71! copy density to the global arrays
72do ias=1,natmtot
73 is=idxis(ias)
74 rhomt(1:npcmt(is),ias)=rhomt_(1:npcmt(is),ias)
75end do
76rhoir(1:ngtc)=rhoir_(1:ngtc)
77if (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
85end if
86deallocate(rhomt_,rhoir_,magmt_,magir_)
87! convert muffin-tin density/magnetisation to spherical harmonics
88call rhomagsh
89call 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
98call rfmtctof(rhomt)
99! convert the interstitial density from coarse to fine grid
100call rfirctof(rhoir,rhoir)
101!$OMP SECTION
102if (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
114end if
115!$OMP END PARALLEL SECTIONS
116call freethd(nthd)
117end subroutine
118
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition getevecfv.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition match.f90:10
real(8), dimension(:,:,:,:), allocatable vgkc
Definition modmain.f90:505
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
integer nfgrzc
Definition modmain.f90:414
integer ngtot
Definition modmain.f90:390
logical ncmag
Definition modmain.f90:240
real(8), dimension(:,:,:), allocatable gkc
Definition modmain.f90:507
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
integer, dimension(3) ngdgc
Definition modmain.f90:388
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
character(256) filext
Definition modmain.f90:1300
logical spinpol
Definition modmain.f90:228
integer, dimension(:,:), allocatable ngk
Definition modmain.f90:497
integer, dimension(:,:,:), allocatable igkig
Definition modmain.f90:501
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer lmmaxapw
Definition modmain.f90:199
integer apwordmax
Definition modmain.f90:760
integer ngtc
Definition modmain.f90:392
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer nspnfv
Definition modmain.f90:289
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
integer npcmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer nmatmax
Definition modmain.f90:848
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), pointer, contiguous magir
Definition modmain.f90:616
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
integer, dimension(:), allocatable igrzfc
Definition modmain.f90:418
integer ngvc
Definition modmain.f90:398
integer nstfv
Definition modmain.f90:884
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition modmain.f90:509
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
integer ndmag
Definition modmain.f90:238
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rfirctof(rfirc, rfir)
Definition rfirctof.f90:10
subroutine rfmtctof(rfmt)
Definition rfmtctof.f90:10
subroutine rhomagk(ngp, igpig, wppt, occsvp, apwalm, evecfv, evecsv, rhomt_, rhoir_, magmt_, magir_)
Definition rhomagk.f90:11
subroutine rhomagsh
Definition rhomagsh.f90:10
subroutine rhomagv
Definition rhomagv.f90:7
subroutine symrf(nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld, rfmt, rfir)
Definition symrf.f90:11
subroutine symrvf(tspin, tnc, nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld1, rvfmt, ld2, rvfir)
Definition symrvf.f90:11