The Elk Code
 
Loading...
Searching...
No Matches
drhomagk.f90
Go to the documentation of this file.
1
2! Copyright (C) 2012 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 drhomagk(ngp,ngpq,igpig,igpqig,occsvp,doccsvp,apwalm,apwalmq, &
7 dapwalm,evecfv,devecfv,evecsv,devecsv)
8use modmain
9use modphonon
10implicit none
11! arguments
12integer, intent(in) :: ngp(nspnfv),ngpq(nspnfv)
13integer, intent(in) :: igpig(ngkmax,nspnfv),igpqig(ngkmax,nspnfv)
14real(8), intent(in) :: occsvp(nstsv),doccsvp(nstsv)
15complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
16complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
17complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw,nspnfv)
18complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv)
19complex(8), intent(in) :: devecfv(nmatmax,nstfv,nspnfv)
20complex(8), intent(in) :: evecsv(nstsv,nstsv),devecsv(nstsv,nstsv)
21! local variables
22integer nst,ist,jst
23integer is,ias,npc
24real(8) wo,dwo
25! automatic arrays
26integer idx(nstsv)
27! allocatable arrays
28complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
29complex(8), allocatable :: dwfmt(:,:,:,:),dwfir(:,:,:)
30! count and index the occupied states
31nst=0
32do ist=1,nstsv
33 if (abs(occsvp(ist)) < epsocc) cycle
34 nst=nst+1
35 idx(nst)=ist
36end do
37! generate the wavefunctions
38allocate(wfmt(npcmtmax,natmtot,nspinor,nst),wfir(ngtc,nspinor,nst))
39call genwfsv(.false.,.false.,nst,idx,ngdgc,igfc,ngp,igpig,apwalm,evecfv,evecsv,&
40 wfmt,ngtc,wfir)
41! generate the wavefunction derivatives
42allocate(dwfmt(npcmtmax,natmtot,nspinor,nst),dwfir(ngtc,nspinor,nst))
43call gendwfsv(.false.,.false.,nst,idx,ngp,ngpq,igpqig,apwalmq,dapwalm,evecfv, &
44 devecfv,evecsv,devecsv,dwfmt,ngtc,dwfir)
45! loop over occupied states
46do ist=1,nst
47 jst=idx(ist)
48 wo=2.d0*wkptnr*occsvp(jst)
49 dwo=wkptnr*doccsvp(jst)
50!----------------------------------------------!
51! muffin-tin density and magnetisation !
52!----------------------------------------------!
53 do ias=1,natmtot
54 is=idxis(ias)
55 npc=npcmt(is)
56!$OMP CRITICAL(drhomagk_1)
57 if (spinpol) then
58! spin-polarised
59 if (ncmag) then
60! non-collinear
61 call drmk1(npc,wo,wfmt(:,ias,1,jst),wfmt(:,ias,2,jst), &
62 dwfmt(:,ias,1,jst),dwfmt(:,ias,2,jst),drhomt(:,ias),dmagmt(:,ias,1), &
63 dmagmt(:,ias,2),dmagmt(:,ias,3))
64 if (tphq0) then
65 call drmk01(npc,dwo,wfmt(:,ias,1,jst),wfmt(:,ias,2,jst), &
66 drhomt(:,ias),dmagmt(:,ias,1),dmagmt(:,ias,2),dmagmt(:,ias,3))
67 end if
68 else
69! collinear
70 call drmk2(npc,wo,wfmt(:,ias,1,jst),wfmt(:,ias,2,jst), &
71 dwfmt(:,ias,1,jst),dwfmt(:,ias,2,jst),drhomt(:,ias),dmagmt(:,ias,1))
72 if (tphq0) then
73 call drmk02(npc,dwo,wfmt(:,ias,1,jst),wfmt(:,ias,2,jst), &
74 drhomt(:,ias),dmagmt(:,ias,1))
75 end if
76 end if
77 else
78! spin-unpolarised
79 call drmk3(npc,wo,wfmt(:,ias,1,jst),dwfmt(:,ias,1,jst),drhomt(:,ias))
80 if (tphq0) then
81 call drmk03(npc,dwo,wfmt(:,ias,1,jst),drhomt(:,ias))
82 end if
83 end if
84!$OMP END CRITICAL(drhomagk_1)
85 end do
86!------------------------------------------------!
87! interstitial density and magnetisation !
88!------------------------------------------------!
89!$OMP CRITICAL(drhomagk_2)
90 if (spinpol) then
91! spin-polarised
92 if (ncmag) then
93 call drmk1(ngtc,wo,wfir(:,1,jst),wfir(:,2,jst),dwfir(:,1,jst), &
94 dwfir(:,2,jst),drhoir,dmagir,dmagir(:,2),dmagir(:,3))
95 if (tphq0) then
96 call drmk01(ngtc,dwo,wfir(:,1,jst),wfir(:,2,jst),drhoir,dmagir, &
97 dmagir(:,2),dmagir(:,3))
98 end if
99 else
100! collinear
101 call drmk2(ngtc,wo,wfir(:,1,jst),wfir(:,2,jst),dwfir(:,1,jst), &
102 dwfir(:,2,jst),drhoir,dmagir)
103 if (tphq0) then
104 call drmk02(ngtc,dwo,wfir(:,1,jst),wfir(:,2,jst),drhoir,dmagir)
105 end if
106 end if
107 else
108! spin-unpolarised
109 call drmk3(ngtc,wo,wfir(:,1,jst),dwfir(:,1,jst),drhoir)
110 if (tphq0) then
111 call drmk03(ngtc,dwo,wfir(:,1,jst),drhoir)
112 end if
113 end if
114!$OMP END CRITICAL(drhomagk_2)
115! end loop over states
116end do
117deallocate(wfmt,wfir,dwfmt,dwfir)
118return
119
120contains
121
122pure subroutine drmk1(n,wo,wf1,wf2,dwf1,dwf2,drho,dmag1,dmag2,dmag3)
123implicit none
124! arguments
125integer, intent(in) :: n
126real(8), intent(in) :: wo
127complex(8), intent(in) :: wf1(n),wf2(n)
128complex(8), intent(in) :: dwf1(n),dwf2(n)
129complex(8), intent(inout) :: drho(n)
130complex(8), intent(inout) :: dmag1(n),dmag2(n),dmag3(n)
131! local variables
132integer i
133complex(8) z1,z2,z3,z4,z5,z6
134do i=1,n
135 z1=conjg(wf1(i))
136 z2=conjg(wf2(i))
137 z3=dwf1(i)
138 z4=dwf2(i)
139 z5=z1*z3
140 z6=z2*z4
141 drho(i)=drho(i)+wo*(z5+z6)
142 dmag3(i)=dmag3(i)+wo*(z5-z6)
143 z5=z1*z4
144 z6=z2*z3
145 dmag1(i)=dmag1(i)+wo*(z5+z6)
146 z5=z5-z6
147 dmag2(i)=dmag2(i)+wo*zmi*z5
148end do
149end subroutine
150
151pure subroutine drmk01(n,dwo,wf1,wf2,drho,dmag1,dmag2,dmag3)
152implicit none
153! arguments
154integer, intent(in) :: n
155real(8), intent(in) :: dwo
156complex(8), intent(in) :: wf1(n),wf2(n)
157complex(8), intent(inout) :: drho(n)
158complex(8), intent(inout) :: dmag1(n),dmag2(n),dmag3(n)
159! local variables
160integer i
161real(8) t1,t2
162complex(8) z1,z2
163do i=1,n
164 z1=wf1(i)
165 z2=wf2(i)
166 t1=z1%re**2+z1%im**2
167 t2=z2%re**2+z2%im**2
168 z1=conjg(z1)*z2
169 drho(i)=drho(i)+dwo*(t1+t2)
170 dmag1(i)=dmag1(i)+dwo*2.d0*z1%re
171 dmag2(i)=dmag2(i)+dwo*2.d0*z1%im
172 dmag3(i)=dmag3(i)+dwo*(t1-t2)
173end do
174end subroutine
175
176pure subroutine drmk2(n,wo,wf1,wf2,dwf1,dwf2,drho,dmag)
177implicit none
178! arguments
179integer, intent(in) :: n
180real(8), intent(in) :: wo
181complex(8), intent(in) :: wf1(n),wf2(n)
182complex(8), intent(in) :: dwf1(n),dwf2(n)
183complex(8), intent(inout) :: drho(n),dmag(n)
184! local variables
185integer i
186complex(8) z1,z2
187do i=1,n
188 z1=conjg(wf1(i))*dwf1(i)
189 z2=conjg(wf2(i))*dwf2(i)
190 drho(i)=drho(i)+wo*(z1+z2)
191 dmag(i)=dmag(i)+wo*(z1-z2)
192end do
193end subroutine
194
195pure subroutine drmk02(n,dwo,wf1,wf2,drho,dmag)
196implicit none
197! arguments
198integer, intent(in) :: n
199real(8), intent(in) :: dwo
200complex(8), intent(in) :: wf1(n),wf2(n)
201complex(8), intent(inout) :: drho(n),dmag(n)
202! local variables
203integer i
204real(8) t1,t2
205do i=1,n
206 t1=dble(wf1(i))**2+aimag(wf1(i))**2
207 t2=dble(wf2(i))**2+aimag(wf2(i))**2
208 drho(i)=drho(i)+dwo*(t1+t2)
209 dmag(i)=dmag(i)+dwo*(t1-t2)
210end do
211end subroutine
212
213pure subroutine drmk3(n,wo,wf,dwf,drho)
214implicit none
215! arguments
216integer, intent(in) :: n
217real(8), intent(in) :: wo
218complex(8), intent(in) :: wf(n),dwf(n)
219complex(8), intent(inout) :: drho(n)
220drho(1:n)=drho(1:n)+wo*conjg(wf(1:n))*dwf(1:n)
221end subroutine
222
223pure subroutine drmk03(n,dwo,wf,drho)
224implicit none
225! arguments
226integer, intent(in) :: n
227real(8), intent(in) :: dwo
228complex(8), intent(in) :: wf(n)
229complex(8), intent(inout) :: drho(n)
230drho(1:n)=drho(1:n)+dwo*(dble(wf(1:n))**2+aimag(wf(1:n))**2)
231end subroutine
232
233end subroutine
234
subroutine drhomagk(ngp, ngpq, igpig, igpqig, occsvp, doccsvp, apwalm, apwalmq, dapwalm, evecfv, devecfv, evecsv, devecsv)
Definition drhomagk.f90:8
pure subroutine drmk3(n, wo, wf, dwf, drho)
Definition drhomagk.f90:214
pure subroutine drmk1(n, wo, wf1, wf2, dwf1, dwf2, drho, dmag1, dmag2, dmag3)
Definition drhomagk.f90:123
pure subroutine drmk2(n, wo, wf1, wf2, dwf1, dwf2, drho, dmag)
Definition drhomagk.f90:177
pure subroutine drmk01(n, dwo, wf1, wf2, drho, dmag1, dmag2, dmag3)
Definition drhomagk.f90:152
pure subroutine drmk02(n, dwo, wf1, wf2, drho, dmag)
Definition drhomagk.f90:196
pure subroutine drmk03(n, dwo, wf, drho)
Definition drhomagk.f90:224
subroutine gendwfsv(tsh, tgp, nst, idx, ngp, ngpq, igpqig, apwalmq, dapwalm, evecfv, devecfv, evecsv, devecsv, dwfmt, ld, dwfir)
Definition gendwfsv.f90:8
subroutine genwfsv(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition genwfsv.f90:11
logical ncmag
Definition modmain.f90:240
complex(8), parameter zmi
Definition modmain.f90:1239
integer nspinor
Definition modmain.f90:267
integer, dimension(3) ngdgc
Definition modmain.f90:388
logical spinpol
Definition modmain.f90:228
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer ngtc
Definition modmain.f90:392
integer npcmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
real(8) epsocc
Definition modmain.f90:900
real(8) wkptnr
Definition modmain.f90:477
logical tphq0
Definition modphonon.f90:17
complex(8), dimension(:,:,:), allocatable dmagmt
complex(8), dimension(:,:), allocatable dmagir
complex(8), dimension(:), allocatable drhoir
Definition modphonon.f90:98
complex(8), dimension(:,:), allocatable drhomt
Definition modphonon.f90:98