The Elk Code
rhomaguk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2017 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 
6 subroutine rhomaguk(ik0,lock,evecu)
7 use modmain
8 use modulr
9 use modomp
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik0
13 integer(omp_lock_kind), intent(inout) :: lock(nqpt)
14 complex(8), intent(in) :: evecu(nstulr,nstulr)
15 ! local variables
16 integer ik,ikpa,jkpa
17 integer nst,ist,jst,i,j
18 integer ngk0,is,ias
19 integer npc,ir,nthd
20 real(8) ts0,ts1
21 real(4) wo
22 ! automatic arrays
23 integer idx(nstsv)
24 complex(8) zfft(nqpt)
25 ! allocatable arrays
26 complex(8), allocatable :: apwalm(:,:,:,:),evecfv(:,:),evecsv(:,:)
27 complex(8), allocatable :: evectv(:,:,:),evecsvt(:,:)
28 complex(4), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
29 call timesec(ts0)
30 ! central k-point
31 ik=(ik0-1)*nkpa+1
32 ! number of G+k-vectors for central k-point
33 ngk0=ngk(1,ik)
34 ! get the eigenvectors from file
35 allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
36 call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
37 call getevecsv(filext,ik,vkl(:,ik),evecsv)
38 ! find the matching coefficients
39 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot))
40 call match(ngk0,vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),apwalm)
41 allocate(evectv(nstsv,nstsv,nqpt))
42 call holdthd(nqpt,nthd)
43 !$OMP PARALLEL DEFAULT(SHARED) &
44 !$OMP PRIVATE(zfft,evecsvt,wfmt,wfir) &
45 !$OMP PRIVATE(ikpa,jkpa,ist,jst,i,j) &
46 !$OMP PRIVATE(ir,wo,ias,is,npc) &
47 !$OMP NUM_THREADS(nthd)
48 allocate(evecsvt(nstsv,nstsv))
49 allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfir(ngtc,nspinor,nstsv))
50 ! loop over long-range states in subsets of size nstsv
51 do jkpa=1,nkpa
52 ! number of and index to occupied states in subset
53 !$OMP BARRIER
54 !$OMP SINGLE
55  nst=0
56  do jst=1,nstsv
57  j=(jkpa-1)*nstsv+jst
58  if (abs(occulr(j,ik0)) < epsocc) cycle
59  nst=nst+1
60  idx(nst)=jst
61  end do
62 !$OMP END SINGLE
63  if (nst == 0) cycle
64 !$OMP DO SCHEDULE(DYNAMIC)
65  do jst=1,nst
66  j=(jkpa-1)*nstsv+idx(jst)
67  do ist=1,nstsv
68  zfft(1:nqpt)=0.d0
69  do ikpa=1,nkpa
70  i=(ikpa-1)*nstsv+ist
71 ! store the long-range state in FFT Q-space
72  zfft(iqfft(ikpa))=evecu(i,j)
73  end do
74 ! Fourier transform to R-space
75  call zfftifc(3,ngridq,1,zfft)
76  evectv(ist,jst,1:nqpt)=zfft(1:nqpt)
77  end do
78  end do
79 !$OMP END DO
80 ! parallel loop over R-points
81 !$OMP DO SCHEDULE(DYNAMIC)
82  do ir=1,nqpt
83 ! convert third-variational states to second-variational states
84  call zgemm('N','N',nstsv,nst,nstsv,zone,evecsv,nstsv,evectv(:,:,ir), &
85  nstsv,zzero,evecsvt,nstsv)
86 ! generate the wavefunctions in single-precision
87  call genwfsv_sp(.false.,.false.,nst,[0],ngdgc,igfc,ngk0,igkig(:,1,ik), &
88  apwalm,evecfv,evecsvt,wfmt,ngtc,wfir)
89 ! loop over second-variational states
90  do jst=1,nst
91  j=(jkpa-1)*nstsv+idx(jst)
92  wo=occulr(j,ik0)*wkpt(ik)
93 ! add to the density and magnetisation
94  call omp_set_lock(lock(ir))
95 ! muffin-tin part
96  do ias=1,natmtot
97  is=idxis(ias)
98  npc=npcmt(is)
99  if (spinpol) then
100  if (ncmag) then
101  call rmk1(npc,wo,wfmt(:,ias,1,jst),wfmt(:,ias,2,jst), &
102  rhormt(:,ias,ir),magrmt(:,ias,1,ir),magrmt(:,ias,2,ir), &
103  magrmt(:,ias,3,ir))
104  else
105  call rmk2(npc,wo,wfmt(:,ias,1,jst),wfmt(:,ias,2,jst), &
106  rhormt(:,ias,ir),magrmt(:,ias,1,ir))
107  end if
108  else
109  call rmk3(npc,wo,wfmt(:,ias,1,jst),rhormt(:,ias,ir))
110  end if
111  end do
112 ! interstitial part
113  if (spinpol) then
114  if (ncmag) then
115  call rmk1(ngtc,wo,wfir(:,1,jst),wfir(:,2,jst),rhorir(:,ir), &
116  magrir(:,1,ir),magrir(:,2,ir),magrir(:,3,ir))
117  else
118  call rmk2(ngtc,wo,wfir(:,1,jst),wfir(:,2,jst),rhorir(:,ir), &
119  magrir(:,1,ir))
120  end if
121  else
122  call rmk3(ngtc,wo,wfir(:,1,jst),rhorir(:,ir))
123  end if
124  call omp_unset_lock(lock(ir))
125  end do
126 ! end parallel loop over R-points
127  end do
128 !$OMP END DO
129 end do
130 deallocate(evecsvt,wfmt,wfir)
131 !$OMP END PARALLEL
132 call freethd(nthd)
133 deallocate(apwalm,evecfv,evecsv,evectv)
134 call timesec(ts1)
135 !$OMP ATOMIC
136 timerho=timerho+ts1-ts0
137 
138 contains
139 
140 pure subroutine rmk1(n,wo,wf1,wf2,rho,mag1,mag2,mag3)
141 implicit none
142 ! arguments
143 integer, intent(in) :: n
144 real(4), intent(in) :: wo
145 complex(4), intent(in) :: wf1(n),wf2(n)
146 real(8), intent(inout) :: rho(n),mag1(n),mag2(n),mag3(n)
147 ! local variables
148 integer i
149 real(4) wo2,a1,b1,a2,b2,t1,t2
150 wo2=2.e0*wo
151 !$OMP SIMD PRIVATE(a1,b1,a2,b2,t1,t2) SIMDLEN(8)
152 do i=1,n
153  a1=real(wf1(i)); b1=aimag(wf1(i))
154  a2=real(wf2(i)); b2=aimag(wf2(i))
155  t1=a1**2+b1**2; t2=a2**2+b2**2
156  mag1(i)=mag1(i)+wo2*(a1*a2+b1*b2)
157  mag2(i)=mag2(i)+wo2*(a1*b2-b1*a2)
158  mag3(i)=mag3(i)+wo*(t1-t2)
159  rho(i)=rho(i)+wo*(t1+t2)
160 end do
161 end subroutine
162 
163 pure subroutine rmk2(n,wo,wf1,wf2,rho,mag)
164 implicit none
165 ! arguments
166 integer, intent(in) :: n
167 real(4), intent(in) :: wo
168 complex(4), intent(in) :: wf1(n),wf2(n)
169 real(8), intent(inout) :: rho(n),mag(n)
170 ! local variables
171 integer i
172 real(4) t1,t2
173 !$OMP SIMD PRIVATE(t1,t2) SIMDLEN(8)
174 do i=1,n
175  t1=real(wf1(i))**2+aimag(wf1(i))**2
176  t2=real(wf2(i))**2+aimag(wf2(i))**2
177  mag(i)=mag(i)+wo*(t1-t2)
178  rho(i)=rho(i)+wo*(t1+t2)
179 end do
180 end subroutine
181 
182 pure subroutine rmk3(n,wo,wf,rho)
183 implicit none
184 ! arguments
185 integer, intent(in) :: n
186 real(4), intent(in) :: wo
187 complex(4), intent(in) :: wf(n)
188 real(8), intent(inout) :: rho(n)
189 rho(1:n)=rho(1:n)+wo*(real(wf(1:n))**2+aimag(wf(1:n))**2)
190 end subroutine
191 
192 end subroutine
193 
integer nmatmax
Definition: modmain.f90:853
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
integer npcmtmax
Definition: modmain.f90:216
character(256) filext
Definition: modmain.f90:1301
integer ngtc
Definition: modmain.f90:392
logical spinpol
Definition: modmain.f90:228
integer lmmaxapw
Definition: modmain.f90:199
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) epsocc
Definition: modmain.f90:903
complex(8), parameter zone
Definition: modmain.f90:1240
real(8) timerho
Definition: modmain.f90:1223
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
integer, dimension(:), allocatable iqfft
Definition: modmain.f90:537
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
subroutine zfftifc(nd, n, sgn, z)
Definition: zfftifc_fftw.f90:7
real(8), dimension(:), allocatable wkpt
Definition: modmain.f90:475
real(8), dimension(:,:,:,:), pointer, contiguous magrmt
Definition: modulr.f90:54
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv_sp.f90:8
real(8), dimension(:,:), allocatable occulr
Definition: modulr.f90:99
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:), pointer, contiguous rhormt
Definition: modulr.f90:53
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
real(8), dimension(:,:), pointer, contiguous rhorir
Definition: modulr.f90:53
complex(8), parameter zzero
Definition: modmain.f90:1240
integer, dimension(3) ngridq
Definition: modmain.f90:515
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer apwordmax
Definition: modmain.f90:760
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
pure subroutine rmk3(n, wo, wf, rho)
Definition: rhomagk.f90:215
subroutine timesec(ts)
Definition: timesec.f90:10
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer, dimension(3) ngdgc
Definition: modmain.f90:388
subroutine rhomaguk(ik0, lock, evecu)
Definition: rhomaguk.f90:7
real(8), dimension(:,:,:), pointer, contiguous magrir
Definition: modulr.f90:54
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer natmtot
Definition: modmain.f90:40
pure subroutine rmk1(n, wo, wf1, wf2, rho, mag1, mag2, mag3)
Definition: rhomagk.f90:172
pure subroutine rmk2(n, wo, wf1, wf2, rho, mag)
Definition: rhomagk.f90:196
Definition: modulr.f90:6
logical ncmag
Definition: modmain.f90:240
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nkpa
Definition: modulr.f90:24
integer nstfv
Definition: modmain.f90:887