The Elk Code
 
Loading...
Searching...
No Matches
potxcu.f90
Go to the documentation of this file.
1
2! Copyright (C) 2018 T. Mueller, 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 potxcu
7use modmain
8use modulr
9use modmpi
10use modomp
11implicit none
12! local variables
13integer ifq,idm,is,ias
14integer ir,npc,n,lp,nthd
15complex(8) z1,z2
16! allocatable arrays
17real(8), allocatable :: vxcrmt(:,:,:),vxcrir(:,:)
18real(8), allocatable :: bxcrmt(:,:,:,:),bxcrir(:,:,:)
19real(8), allocatable :: rhomt_(:,:),rhoir_(:)
20real(8), allocatable :: magmt_(:,:,:),magir_(:,:)
21real(8), allocatable :: vxcmt_(:,:),bxcmt_(:,:,:)
22complex(8), allocatable :: vxcqmt(:,:,:),vxcqir(:,:)
23allocate(vxcrmt(npcmtmax,natmtot,nqpt),vxcrir(ngtot,nqpt))
24if (spinpol) then
25 allocate(bxcrmt(npcmtmax,natmtot,ndmag,nqpt))
26 allocate(bxcrir(ngtot,ndmag,nqpt))
27end if
28call holdthd(nqpt/np_mpi,nthd)
29!$OMP PARALLEL DEFAULT(SHARED) &
30!$OMP PRIVATE(rhomt_,rhoir_,vxcmt_) &
31!$OMP PRIVATE(magmt_,magir_,bxcmt_) &
32!$OMP PRIVATE(ias,is,idm) &
33!$OMP NUM_THREADS(nthd)
34allocate(rhomt_(npmtmax,natmtot),rhoir_(ngtot))
35allocate(vxcmt_(npmtmax,natmtot))
36if (spinpol) then
37 allocate(magmt_(npmtmax,natmtot,ndmag),magir_(ngtot,ndmag))
38 allocate(bxcmt_(npmtmax,natmtot,ndmag))
39end if
40!$OMP DO SCHEDULE(DYNAMIC)
41do ir=1,nqpt
42! distribute among MPI processes
43 if (mod(ir-1,np_mpi) /= lp_mpi) cycle
44! convert muffin-tin density and magnetisation from coarse to fine radial mesh
45 do ias=1,natmtot
46 is=idxis(ias)
47 rhomt_(1:npcmt(is),ias)=rhormt(1:npcmt(is),ias,ir)
48 end do
49 call rfmtctof(rhomt_)
50 do idm=1,ndmag
51 do ias=1,natmtot
52 is=idxis(ias)
53 magmt_(1:npcmt(is),ias,idm)=magrmt(1:npcmt(is),ias,idm,ir)
54 end do
55 call rfmtctof(magmt_(:,:,idm))
56 end do
57! convert interstitial density and magnetisation from coarse to fine grid
58 call rfirctof(rhorir(:,ir),rhoir_)
59 do idm=1,ndmag
60 call rfirctof(magrir(:,idm,ir),magir_(:,idm))
61 end do
62! calculate the exchange-correlation potential and magnetic field
63 call potxc(.false.,xctype,rhomt_,rhoir_,magmt_,magir_,taumt,tauir,exmt,exir, &
64 ecmt,ecir,vxcmt_,vxcrir(:,ir),bxcmt_,bxcrir(:,:,ir),wxcmt,wxcir)
65! convert muffin-tin potential and field from fine to coarse radial mesh
66 do ias=1,natmtot
67 is=idxis(ias)
68 call rfmtftoc(nrcmt(is),nrcmti(is),vxcmt_(:,ias),vxcrmt(:,ias,ir))
69 end do
70 do idm=1,ndmag
71 do ias=1,natmtot
72 is=idxis(ias)
73 call rfmtftoc(nrcmt(is),nrcmti(is),bxcmt_(:,ias,idm),bxcrmt(:,ias,idm,ir))
74 end do
75 end do
76end do
77!$OMP END DO
78deallocate(rhomt_,rhoir_,vxcmt_)
79if (spinpol) deallocate(magmt_,magir_,bxcmt_)
80!$OMP END PARALLEL
81call freethd(nthd)
82! broadcast potentials and fields to every MPI process
83if (np_mpi > 1) then
85 do ir=1,nqpt
86 lp=mod(ir-1,np_mpi)
87 call mpi_bcast(vxcrmt(:,:,ir),n,mpi_double_precision,lp,mpicom,ierror)
88 end do
89 do ir=1,nqpt
90 lp=mod(ir-1,np_mpi)
91 call mpi_bcast(vxcrir(:,ir),ngtot,mpi_double_precision,lp,mpicom,ierror)
92 end do
93 if (spinpol) then
95 do ir=1,nqpt
96 lp=mod(ir-1,np_mpi)
97 call mpi_bcast(bxcrmt(:,:,:,ir),n,mpi_double_precision,lp,mpicom,ierror)
98 end do
99 n=ngtot*ndmag
100 do ir=1,nqpt
101 lp=mod(ir-1,np_mpi)
102 call mpi_bcast(bxcrir(:,:,ir),n,mpi_double_precision,lp,mpicom,ierror)
103 end do
104 end if
105end if
106allocate(vxcqmt(npcmtmax,natmtot,nfqrz),vxcqir(ngtot,nfqrz))
107! Fourier transform exchange-correlation potential to Q-space
108call rfzfftq(-1,1,ngtot,vxcrmt,vxcrir,vxcqmt,vxcqir)
109deallocate(vxcrmt,vxcrir)
110! add V_xc and external Coulomb potential to Kohn-Sham potential
111call holdthd(2*nfqrz,nthd)
112!$OMP PARALLEL DEFAULT(SHARED) &
113!$OMP PRIVATE(z1,ias,is,npc) &
114!$OMP NUM_THREADS(nthd)
115!$OMP DO SCHEDULE(DYNAMIC)
116do ifq=1,nfqrz
117 z1=vclq(ifq)
118 do ias=1,natmtot
119 is=idxis(ias)
120 npc=npcmt(is)
121 vsqmt(1:npc,ias,ifq)=vsqmt(1:npc,ias,ifq)+vxcqmt(1:npc,ias,ifq)+z1
122 end do
123end do
124!$OMP END DO NOWAIT
125!$OMP DO SCHEDULE(DYNAMIC)
126do ifq=1,nfqrz
127 z1=vclq(ifq)
128 vsqir(1:ngtot,ifq)=vsqir(1:ngtot,ifq)+vxcqir(1:ngtot,ifq)+z1
129end do
130!$OMP END DO
131!$OMP END PARALLEL
132call freethd(nthd)
133deallocate(vxcqmt,vxcqir)
134if (spinpol) then
135! Fourier transform the exchange-correlation magnetic field to Q-space
136 call rfzfftq(-1,ndmag,ngtot,bxcrmt,bxcrir,bsqmt,bsqir)
137 deallocate(bxcrmt,bxcrir)
138! add external magnetic fields
139 call holdthd(2*nfqrz,nthd)
140!$OMP PARALLEL DEFAULT(SHARED) &
141!$OMP PRIVATE(idm,z1,z2,ias,is,npc) &
142!$OMP NUM_THREADS(nthd)
143!$OMP DO SCHEDULE(DYNAMIC)
144 do ifq=1,nfqrz
145 do idm=1,ndmag
146 z1=bfcq(idm,ifq)
147 if (tbdipu) z1=z1+bdipq(idm,ifq)
148 do ias=1,natmtot
149 is=idxis(ias)
150 npc=npcmt(is)
151 z2=z1+bfcmtq(ias,idm,ifq)
152 bsqmt(1:npc,ias,idm,ifq)=bsqmt(1:npc,ias,idm,ifq)+z2
153 end do
154 end do
155 end do
156!$OMP END DO NOWAIT
157!$OMP DO SCHEDULE(DYNAMIC)
158 do ifq=1,nfqrz
159 do idm=1,ndmag
160 z1=bfcq(idm,ifq)
161 if (tbdipu) z1=z1+bdipq(idm,ifq)
162 bsqir(1:ngtot,idm,ifq)=bsqir(1:ngtot,idm,ifq)+z1
163 end do
164 end do
165!$OMP END DO
166!$OMP END PARALLEL
167 call freethd(nthd)
168end if
169end subroutine
170
integer ngtot
Definition modmain.f90:390
integer nfqrz
Definition modmain.f90:539
real(8), dimension(:,:), allocatable exmt
Definition modmain.f90:630
real(8), dimension(:,:,:), allocatable taumt
Definition modmain.f90:672
logical spinpol
Definition modmain.f90:228
real(8), dimension(:,:), allocatable ecmt
Definition modmain.f90:632
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(3) xctype
Definition modmain.f90:588
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
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 npmtmax
Definition modmain.f90:216
real(8), dimension(:), allocatable exir
Definition modmain.f90:630
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
real(8), dimension(:,:), allocatable tauir
Definition modmain.f90:672
real(8), dimension(:), allocatable wxcir
Definition modmain.f90:676
integer ndmag
Definition modmain.f90:238
real(8), dimension(:,:), allocatable wxcmt
Definition modmain.f90:676
real(8), dimension(:), allocatable ecir
Definition modmain.f90:632
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
complex(8), dimension(:), allocatable vclq
Definition modulr.f90:65
logical tbdipu
Definition modulr.f90:78
complex(8), dimension(:,:,:,:), pointer, contiguous bsqmt
Definition modulr.f90:84
complex(8), dimension(:,:), allocatable bfcq
Definition modulr.f90:70
real(8), dimension(:,:,:), allocatable rhormt
Definition modulr.f90:52
real(8), dimension(:,:,:,:), allocatable magrmt
Definition modulr.f90:53
complex(8), dimension(:,:,:), pointer, contiguous vsqmt
Definition modulr.f90:83
complex(8), dimension(:,:,:), pointer, contiguous bsqir
Definition modulr.f90:84
complex(8), dimension(:,:), allocatable bdipq
Definition modulr.f90:80
real(8), dimension(:,:,:), allocatable magrir
Definition modulr.f90:53
complex(8), dimension(:,:), pointer, contiguous vsqir
Definition modulr.f90:83
real(8), dimension(:,:), allocatable rhorir
Definition modulr.f90:52
complex(8), dimension(:,:,:), allocatable bfcmtq
Definition modulr.f90:72
subroutine potxc(tsh, xctype_, rhomt_, rhoir_, magmt_, magir_, taumt_, tauir_, exmt_, exir_, ecmt_, ecir_, vxcmt_, vxcir_, bxcmt_, bxcir_, wxcmt_, wxcir_)
Definition potxc.f90:11
subroutine potxcu
Definition potxcu.f90:7
subroutine rfirctof(rfirc, rfir)
Definition rfirctof.f90:10
subroutine rfmtctof(rfmt)
Definition rfmtctof.f90:10
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition rfmtftoc.f90:7
subroutine rfzfftq(sgn, nf, ngt, rfmt, rfir, zfmt, zfir)
Definition rfzfftq.f90:7