The Elk Code
rhomagk.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2010 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: rhomagk
8 ! !INTERFACE:
9 subroutine rhomagk(ngp,igpig,wppt,occsvp,apwalm,evecfv,evecsv,rhomt_,rhoir_, &
10  magmt_,magir_)
11 ! !USES:
12 use modmain
13 use modomp
14 ! !INPUT/OUTPUT PARAMETERS:
15 ! ngp : number of G+p-vectors (in,integer(nspnfv))
16 ! igpig : index from G+p-vectors to G-vectors (in,integer(ngkmax,nspnfv))
17 ! wppt : weight of input p-point (in,real)
18 ! occsvp : occupation number for each state (in,real(nstsv))
19 ! apwalm : APW matching coefficients
20 ! (in,complex(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
21 ! evecfv : first-variational eigenvectors (in,complex(nmatmax,nstfv,nspnfv))
22 ! evecsv : second-variational eigenvectors (in,complex(nstsv,nstsv))
23 ! rhomt_ : muffin-tin density (inout,real(npcmtmax,natmtot))
24 ! rhoir_ : interstitial density (inout,real(ngtc))
25 ! magmt_ : muffin-tin magnetisation (inout,real(npcmtmax,natmtot,ndmag))
26 ! magir_ : interstitial magnetisation (inout,real(ngtc,ndmag))
27 ! !DESCRIPTION:
28 ! Generates the partial valence charge density and magnetisation from the
29 ! eigenvectors at a particular $k$-point. In the muffin-tin region, the
30 ! wavefunction is obtained in terms of its $(l,m)$-components from both the
31 ! APW and local-orbital functions. Using a backward spherical harmonic
32 ! transform (SHT), the wavefunction is converted to real-space and the density
33 ! obtained from its modulus squared. A similar process is used for the
34 ! interstitial density in which the wavefunction in real-space is obtained
35 ! from a Fourier transform of the APW functions. See routines {\tt wfmtsv},
36 ! {\tt genshtmat} and {\tt eveqn}.
37 !
38 ! !REVISION HISTORY:
39 ! Created April 2003 (JKD)
40 ! Removed conversion to spherical harmonics, January 2009 (JKD)
41 ! Partially de-phased the muffin-tin magnetisation for spin-spirals,
42 ! February 2009 (FC, FB & LN)
43 ! Optimisations, July 2010 (JKD)
44 !EOP
45 !BOC
46 implicit none
47 ! arguments
48 integer, intent(in) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
49 real(8), intent(in) :: wppt,occsvp(nstsv)
50 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
51 complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
52 real(8), intent(inout) :: rhomt_(npcmtmax,natmtot),rhoir_(ngtc)
53 real(8), intent(inout) :: magmt_(npcmtmax,natmtot,ndmag),magir_(ngtc,ndmag)
54 ! local variables
55 integer ispn,jspn,nst,ist
56 integer is,ias,npc,j,k
57 integer n,igp,nthd
58 real(8) wo,ts0,ts1
59 complex(8) z1
60 ! automatic arrays
61 integer idx(nstsv)
62 complex(8) wfir(ngtc,nspinor),wfgp(ngkmax)
63 ! allocatable arrays
64 complex(8), allocatable :: wfmt(:,:,:)
65 call timesec(ts0)
66 ! number of and index to occupied states
67 nst=0
68 do ist=1,nstsv
69  if (abs(occsvp(ist)) < epsocc) cycle
70  nst=nst+1
71  idx(nst)=ist
72 end do
73 !----------------------------------------------!
74 ! muffin-tin density and magnetisation !
75 !----------------------------------------------!
76 call holdthd(max(natmtot,nst),nthd)
77 !$OMP PARALLEL DEFAULT(SHARED) &
78 !$OMP PRIVATE(wfmt,wfir,wfgp) &
79 !$OMP PRIVATE(ias,is,npc,j,k,wo) &
80 !$OMP PRIVATE(ispn,jspn,n,ist,z1,igp) &
81 !$OMP NUM_THREADS(nthd)
82 allocate(wfmt(npcmtmax,nspinor,nst))
83 !$OMP DO SCHEDULE(DYNAMIC)
84 do ias=1,natmtot
85  is=idxis(ias)
86  npc=npcmt(is)
87  call wfmtsv(.false.,lradstp,is,ias,nst,idx,ngp,apwalm,evecfv,evecsv,npcmtmax,&
88  wfmt)
89  do j=1,nst
90  k=idx(j)
91  wo=occsvp(k)*wppt
92 ! add to density and magnetisation
93  if (spinpol) then
94 ! spin-polarised
95  if (ncmag) then
96 ! non-collinear
97  call rmk1(npc,wo,wfmt(:,1,j),wfmt(:,2,j),rhomt_(:,ias),magmt_(:,ias,1),&
98  magmt_(:,ias,2),magmt_(:,ias,3))
99  else
100 ! collinear
101  call rmk2(npc,wo,wfmt(:,1,j),wfmt(:,2,j),rhomt_(:,ias),magmt_(:,ias,1))
102  end if
103  else
104 ! spin-unpolarised
105  call rmk3(npc,wo,wfmt(:,1,j),rhomt_(:,ias))
106  end if
107  end do
108 end do
109 !$OMP END DO NOWAIT
110 deallocate(wfmt)
111 !------------------------------------------------!
112 ! interstitial density and magnetisation !
113 !------------------------------------------------!
114 !$OMP DO ORDERED SCHEDULE(DYNAMIC)
115 do j=1,nst
116  k=idx(j)
117  wo=occsvp(k)*wppt/omega
118  if (tevecsv) then
119 ! generate spinor wavefunction from second-variational eigenvectors
120  do ispn=1,nspinor
121  jspn=jspnfv(ispn)
122  n=ngp(jspn)
123  wfgp(1:n)=0.d0
124  do ist=1,nstfv
125  z1=evecsv((ispn-1)*nstfv+ist,k)
126  if (abs(z1%re)+abs(z1%im) > epswf) then
127  wfgp(1:n)=wfgp(1:n)+z1*evecfv(1:n,ist,jspn)
128  end if
129  end do
130  wfir(1:ngtc,ispn)=0.d0
131  do igp=1,n
132  wfir(igfc(igpig(igp,jspn)),ispn)=wfgp(igp)
133  end do
134 ! Fourier transform wavefunction to real-space
135  call zfftifc(3,ngdgc,1,wfir(:,ispn))
136  end do
137  else
138 ! spin-unpolarised wavefunction
139  wfir(1:ngtc,1)=0.d0
140  do igp=1,ngp(1)
141  wfir(igfc(igpig(igp,1)),1)=evecfv(igp,k,1)
142  end do
143  call zfftifc(3,ngdgc,1,wfir)
144  end if
145 ! add to density and magnetisation
146 !$OMP ORDERED
147  if (spinpol) then
148 ! spin-polarised
149  if (ncmag) then
150 ! non-collinear
151  call rmk1(ngtc,wo,wfir,wfir(:,2),rhoir_,magir_,magir_(:,2),magir_(:,3))
152  else
153 ! collinear
154  call rmk2(ngtc,wo,wfir,wfir(:,2),rhoir_,magir_)
155  end if
156  else
157 ! spin-unpolarised
158  call rmk3(ngtc,wo,wfir,rhoir_)
159  end if
160 !$OMP END ORDERED
161 end do
162 !$OMP END DO
163 !$OMP END PARALLEL
164 call freethd(nthd)
165 call timesec(ts1)
166 !$OMP ATOMIC
167 timerho=timerho+ts1-ts0
168 
169 contains
170 
171 pure subroutine rmk1(n,wo,wf1,wf2,rho,mag1,mag2,mag3)
172 implicit none
173 ! arguments
174 integer, intent(in) :: n
175 real(8), intent(in) :: wo
176 complex(8), intent(in) :: wf1(n),wf2(n)
177 real(8), intent(inout) :: rho(n),mag1(n),mag2(n),mag3(n)
178 ! local variables
179 integer i
180 real(8) wo2,t1,t2
181 real(8) a1,b1,a2,b2
182 wo2=2.d0*wo
183 !$OMP SIMD PRIVATE(a1,b1,a2,b2,t1,t2) SIMDLEN(8)
184 do i=1,n
185  a1=dble(wf1(i)); b1=aimag(wf1(i))
186  a2=dble(wf2(i)); b2=aimag(wf2(i))
187  t1=a1**2+b1**2; t2=a2**2+b2**2
188  mag1(i)=mag1(i)+wo2*(a1*a2+b1*b2)
189  mag2(i)=mag2(i)+wo2*(a1*b2-b1*a2)
190  mag3(i)=mag3(i)+wo*(t1-t2)
191  rho(i)=rho(i)+wo*(t1+t2)
192 end do
193 end subroutine
194 
195 pure subroutine rmk2(n,wo,wf1,wf2,rho,mag)
196 implicit none
197 ! arguments
198 integer, intent(in) :: n
199 real(8), intent(in) :: wo
200 complex(8), intent(in) :: wf1(n),wf2(n)
201 real(8), intent(inout) :: rho(n),mag(n)
202 ! local variables
203 integer i
204 real(8) t1,t2
205 !$OMP SIMD PRIVATE(t1,t2) SIMDLEN(8)
206 do i=1,n
207  t1=dble(wf1(i))**2+aimag(wf1(i))**2
208  t2=dble(wf2(i))**2+aimag(wf2(i))**2
209  mag(i)=mag(i)+wo*(t1-t2)
210  rho(i)=rho(i)+wo*(t1+t2)
211 end do
212 end subroutine
213 
214 pure subroutine rmk3(n,wo,wf,rho)
215 implicit none
216 ! arguments
217 integer, intent(in) :: n
218 real(8), intent(in) :: wo
219 complex(8), intent(in) :: wf(n)
220 real(8), intent(inout) :: rho(n)
221 rho(1:n)=rho(1:n)+wo*(dble(wf(1:n))**2+aimag(wf(1:n))**2)
222 end subroutine
223 
224 end subroutine
225 !EOC
226 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
logical spinpol
Definition: modmain.f90:228
logical tevecsv
Definition: modmain.f90:923
real(8) omega
Definition: modmain.f90:20
Definition: modomp.f90:6
real(8) epsocc
Definition: modmain.f90:903
subroutine rhomagk(ngp, igpig, wppt, occsvp, apwalm, evecfv, evecsv, rhomt_, rhoir_, magmt_, magir_)
Definition: rhomagk.f90:11
real(8) timerho
Definition: modmain.f90:1223
subroutine zfftifc(nd, n, sgn, z)
Definition: zfftifc_fftw.f90:7
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
integer lradstp
Definition: modmain.f90:171
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
pure subroutine rmk3(n, wo, wf, rho)
Definition: rhomagk.f90:215
real(8), parameter epswf
Definition: modmain.f90:845
subroutine timesec(ts)
Definition: timesec.f90:10
integer, dimension(3) ngdgc
Definition: modmain.f90:388
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
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
logical ncmag
Definition: modmain.f90:240
subroutine wfmtsv(tsh, lrstp, is, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, wfmt)
Definition: wfmtsv.f90:7
integer, dimension(2) jspnfv
Definition: modmain.f90:291