The Elk Code
 
Loading...
Searching...
No Matches
gentau.f90
Go to the documentation of this file.
1
2! Copyright (C) 2011 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 gentau
7use modmain
8use modmpi
9use modomp
10implicit none
11! local variables
12integer ik,ispn,is,ias
13integer np,npc,n,nthd
14! allocatable arrays
15real(8), allocatable :: taumt_(:,:,:),tauir_(:,:)
16real(8), allocatable :: rfmt(:,:),rfir(:),rvfmt(:,:,:),rvfir(:,:)
17! tau cannot be computed if wavefunctions do not exist
18if (iscl < 1) then
19 taumt(:,:,:)=0.d0
20 tauir(:,:)=0.d0
21 return
22end if
23! set the kinetic energy density to zero on the coarse grids
24allocate(taumt_(npcmtmax,natmtot,nspinor),tauir_(ngtc,nspinor))
25taumt_(:,:,:)=0.d0
26tauir_(:,:)=0.d0
27call holdthd(nkpt/np_mpi,nthd)
28!$OMP PARALLEL DO DEFAULT(SHARED) &
29!$OMP REDUCTION(+:taumt_,tauir_) &
30!$OMP SCHEDULE(STATIC) &
31!$OMP NUM_THREADS(nthd)
32do ik=1,nkpt
33! distribute among MPI processes
34 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
35 call gentauk(ik,taumt_,tauir_)
36end do
37!$OMP END PARALLEL DO
38call freethd(nthd)
39! add tau from each process and redistribute
40if (np_mpi > 1) then
42 call mpi_allreduce(mpi_in_place,taumt_,n,mpi_double_precision,mpi_sum,mpicom,&
43 ierror)
45 call mpi_allreduce(mpi_in_place,tauir_,n,mpi_double_precision,mpi_sum,mpicom,&
46 ierror)
47end if
48! copy tau to the global arrays
49do ispn=1,nspinor
50 do ias=1,natmtot
51 is=idxis(ias)
52 taumt(1:npcmt(is),ias,ispn)=taumt_(1:npcmt(is),ias,ispn)
53 end do
54end do
55tauir(1:ngtc,1:nspinor)=tauir_(1:ngtc,1:nspinor)
56deallocate(taumt_,tauir_)
57! convert taumt to spherical harmonics
58do ispn=1,nspinor
59 do ias=1,natmtot
60 is=idxis(ias)
61 call rfshtip(nrcmt(is),nrcmti(is),taumt(:,ias,ispn))
62 end do
63end do
64! symmetrise tau
65if (spinpol) then
66! spin-polarised case: convert to scalar-vector form
67 allocate(rfmt(npcmtmax,natmtot),rfir(ngtc))
68 allocate(rvfmt(npcmtmax,natmtot,ndmag),rvfir(ngtc,ndmag))
69 do ias=1,natmtot
70 is=idxis(ias)
71 npc=npcmt(is)
72 rfmt(1:npc,ias)=taumt(1:npc,ias,1)+taumt(1:npc,ias,2)
73 rvfmt(1:npc,ias,1:ndmag-1)=0.d0
74 rvfmt(1:npc,ias,ndmag)=taumt(1:npc,ias,1)-taumt(1:npc,ias,2)
75 end do
76 rfir(1:ngtc)=tauir(1:ngtc,1)+tauir(1:ngtc,2)
77 rvfir(1:ngtc,1:ndmag-1)=0.d0
78 rvfir(1:ngtc,ndmag)=tauir(1:ngtc,1)-tauir(1:ngtc,2)
80 rfmt,rfir)
82 igrzfc,npcmtmax,rvfmt,ngtc,rvfir)
83 do ias=1,natmtot
84 is=idxis(ias)
85 npc=npcmt(is)
86 taumt(1:npc,ias,1)=0.5d0*(rfmt(1:npc,ias)+rvfmt(1:npc,ias,ndmag))
87 taumt(1:npc,ias,2)=0.5d0*(rfmt(1:npc,ias)-rvfmt(1:npc,ias,ndmag))
88 end do
89 tauir(1:ngtc,1)=0.5d0*(rfir(1:ngtc)+rvfir(1:ngtc,ndmag))
90 tauir(1:ngtc,2)=0.5d0*(rfir(1:ngtc)-rvfir(1:ngtc,ndmag))
91 deallocate(rfmt,rfir,rvfmt,rvfir)
92else
93! spin-unpolarised case
96end if
97! convert muffin-tin tau from coarse to fine radial mesh
98do ispn=1,nspinor
99 call rfmtctof(taumt(:,:,ispn))
100end do
101! convert interstitial tau from coarse to fine grid
102do ispn=1,nspinor
103 call rfirctof(tauir(:,ispn),tauir(:,ispn))
104end do
105! generate the core kinetic energy density
106call gentaucr
107do ispn=1,nspinor
108 do ias=1,natmtot
109 is=idxis(ias)
110 np=npmt(is)
111! add the core contribution
112 taumt(1:np,ias,ispn)=taumt(1:np,ias,ispn)+taucr(1:np,ias,ispn)
113! improve stability by smoothing tau
114 call rfmtsm(msmgmt,nrmt(is),nrmti(is),taumt(:,ias,ispn))
115 end do
116end do
117end subroutine
118
subroutine gentau
Definition gentau.f90:7
subroutine gentaucr
Definition gentaucr.f90:7
subroutine gentauk(ik, taumt_, tauir_)
Definition gentauk.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer nfgrzc
Definition modmain.f90:414
logical ncmag
Definition modmain.f90:240
integer nspinor
Definition modmain.f90:267
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer msmgmt
Definition modmain.f90:222
real(8), dimension(:,:,:), allocatable taumt
Definition modmain.f90:672
logical spinpol
Definition modmain.f90:228
real(8), dimension(:,:,:), allocatable taucr
Definition modmain.f90:674
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 ngtc
Definition modmain.f90:392
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer npcmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
integer, dimension(:), allocatable igrzfc
Definition modmain.f90:418
integer ngvc
Definition modmain.f90:398
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer iscl
Definition modmain.f90:1048
real(8), dimension(:,:), allocatable tauir
Definition modmain.f90:672
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 rfmtsm(m, nr, nri, rfmt)
Definition rfmtsm.f90:7
subroutine rfshtip(nr, nri, rfmt)
Definition rfshtip.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