The Elk Code
 
Loading...
Searching...
No Matches
getwfpw.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 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 getwfpw(vpl,vhpl,wfpw)
7use modmain
8use modpw
9implicit none
10! arguments
11real(8), intent(in) :: vpl(3)
12real(8), intent(in) :: vhpl(3,nhkmax,nspnfv)
13complex(8), intent(out) :: wfpw(nhkmax,nspinor,nstsv)
14! local variables
15integer isym,lspl,ilspl,lspn
16integer ik,ist,ihk,ihp,jhp,ig
17integer ispn0,ispn1,jspn,i
18integer recl,nhkmax_,nspinor_,nstsv_
19real(8) vkl_(3),si(3,3)
20real(8) v(3),det,th,t1
21complex(8) su2(2,2),z1,z2
22! automatic arrays
23logical done(nhkmax)
24! allocatable arrays
25complex(8), allocatable :: wfpw_(:,:,:)
26! find the equivalent k-point number and symmetry which rotates vkl to vpl
27call findkpt(vpl,isym,ik)
28! index to spatial rotation in lattice point group
29lspl=lsplsymc(isym)
30! find the record length
31inquire(iolength=recl) vkl_,nhkmax_,nspinor_,nstsv_,wfpw
32!$OMP CRITICAL(u270)
33open(270,file='WFPW.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
34read(270,rec=ik) vkl_,nhkmax_,nspinor_,nstsv_,wfpw
35close(270)
36!$OMP END CRITICAL(u270)
37t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
38if (t1 > epslat) then
39 write(*,*)
40 write(*,'("Error(getwfpw): differing vectors for k-point ",I8)') ik
41 write(*,'(" current : ",3G18.10)') vkl(:,ik)
42 write(*,'(" WFPW.OUT : ",3G18.10)') vkl_
43 write(*,*)
44 stop
45end if
46if (nhkmax /= nhkmax_) then
47 write(*,*)
48 write(*,'("Error(getwfpw): differing nhkmax for k-point ",I8)') ik
49 write(*,'(" current : ",I8)') nhkmax
50 write(*,'(" WFPW.OUT : ",I8)') nhkmax_
51 write(*,*)
52 stop
53end if
54if (nspinor /= nspinor_) then
55 write(*,*)
56 write(*,'("Error(getwfpw): differing nspinor for k-point ",I8)') ik
57 write(*,'(" current : ",I8)') nspinor
58 write(*,'(" WFPW.OUT : ",I8)') nspinor_
59 write(*,*)
60 stop
61end if
62if (nstsv /= nstsv_) then
63 write(*,*)
64 write(*,'("Error(getwfpw): differing nstsv for k-point ",I8)') ik
65 write(*,'(" current : ",I8)') nstsv
66 write(*,'(" WFPW.OUT : ",I8)') nstsv_
67 write(*,*)
68 stop
69end if
70! if p = k then return
71t1=abs(vpl(1)-vkl(1,ik))+abs(vpl(2)-vkl(2,ik))+abs(vpl(3)-vkl(3,ik))
72if (t1 < epslat) return
73!--------------------------------------------------------!
74! translate and rotate wavefunction coefficients !
75!--------------------------------------------------------!
76! allocate temporary copy of wavefunction
77allocate(wfpw_(nhkmax,nspinor,nstsv))
78! the inverse of the spatial symmetry
79ilspl=isymlat(lspl)
80si(:,:)=dble(symlat(:,:,ilspl))
81! loop over first-variational spins
82do jspn=1,nspnfv
83 if (spinsprl) then
84 ispn0=jspn; ispn1=jspn
85 else
86 ispn0=1; ispn1=nspinor
87 end if
88! apply translation operation if required
89 if (tv0symc(isym)) then
90! translation vector is zero
91 do ihk=1,nhk(jspn,ik)
92 wfpw_(ihk,ispn0:ispn1,:)=wfpw(ihk,ispn0:ispn1,:)
93 end do
94 else
95! non-zero translation vector gives a phase factor
96 v(:)=vtcsymc(:,isym)
97 do ihk=1,nhk(jspn,ik)
98 ig=ihkig(ihk,jspn,ik)
99 t1=vgc(1,ig)*v(1)+vgc(2,ig)*v(2)+vgc(3,ig)*v(3)
100 z1=cmplx(cos(t1),-sin(t1),8)
101 wfpw_(ihk,ispn0:ispn1,:)=z1*wfpw(ihk,ispn0:ispn1,:)
102 end do
103 end if
104! apply spatial rotation operation (passive transformation)
105 done(1:nhk(jspn,ik))=.false.
106 i=1
107 do ihk=1,nhk(jspn,ik)
108 call r3mtv(si,vhkl(:,ihk,jspn,ik),v)
109 do ihp=i,nhk(jspn,ik)
110 if (done(ihp)) cycle
111 t1=abs(v(1)-vhpl(1,ihp,jspn)) &
112 +abs(v(2)-vhpl(2,ihp,jspn)) &
113 +abs(v(3)-vhpl(3,ihp,jspn))
114 if (t1 < epslat) then
115 wfpw(ihp,ispn0:ispn1,:)=wfpw_(ihk,ispn0:ispn1,:)
116 done(ihp)=.true.
117 exit
118 end if
119 end do
120 do ihp=i,nhk(jspn,ik)
121 if (.not.done(ihp)) then
122 i=ihp
123 exit
124 end if
125 end do
126 end do
127end do
128! apply spin rotation if required
129if (spinpol) then
130! index to global spin rotation in lattice point group
131 lspn=lspnsymc(isym)
132! if symmetry element is the identity return
133 if (lspn == 1) return
134! find the SU(2) representation of the spin rotation matrix
135 call rotaxang(epslat,symlatc(:,:,lspn),det,v,th)
136 call axangsu2(v,th,su2)
137! apply SU(2) matrix to spinor wavefunctions (active transformation)
138 if (spinsprl) then
139! spin-spiral case
140 wfpw(:,2,:)=0.d0
141 i=1
142 do ihp=1,nhk(1,ik)
143 v(:)=vhpl(:,ihp,1)-vqlss(:)
144 do jhp=i,nhk(2,ik)
145 t1=abs(v(1)-vhpl(1,jhp,2)) &
146 +abs(v(2)-vhpl(2,jhp,2)) &
147 +abs(v(3)-vhpl(3,jhp,2))
148 if (t1 < epslat) then
149 do ist=1,nstsv
150 z1=wfpw(ihp,1,ist)
151 z2=wfpw(jhp,2,ist)
152 wfpw(ihp,1,ist)=su2(1,1)*z1+su2(1,2)*z2
153 wfpw(jhp,2,ist)=su2(2,1)*z1+su2(2,2)*z2
154 end do
155 if (jhp == i) i=i+1
156 goto 10
157 end if
158 end do
159 wfpw(ihp,1,:)=0.d0
16010 continue
161 end do
162 else
163! normal spin case
164 do ist=1,nstsv
165 do ihp=1,nhk(1,ik)
166 z1=wfpw(ihp,1,ist)
167 z2=wfpw(ihp,2,ist)
168 wfpw(ihp,1,ist)=su2(1,1)*z1+su2(1,2)*z2
169 wfpw(ihp,2,ist)=su2(2,1)*z1+su2(2,2)*z2
170 end do
171 end do
172 end if
173end if
174deallocate(wfpw_)
175end subroutine
176
pure subroutine axangsu2(v, th, su2)
Definition axangsu2.f90:10
subroutine findkpt(vpl, isym, ik)
Definition findkpt.f90:7
subroutine getwfpw(vpl, vhpl, wfpw)
Definition getwfpw.f90:7
integer, dimension(48) isymlat
Definition modmain.f90:348
real(8), dimension(3, 3, 48) symlatc
Definition modmain.f90:350
logical spinpol
Definition modmain.f90:228
integer, dimension(maxsymcrys) lspnsymc
Definition modmain.f90:366
real(8) epslat
Definition modmain.f90:24
logical spinsprl
Definition modmain.f90:283
logical, dimension(maxsymcrys) tv0symc
Definition modmain.f90:362
real(8), dimension(3, maxsymcrys) vtcsymc
Definition modmain.f90:360
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
real(8), dimension(3) vqlss
Definition modmain.f90:293
integer, dimension(3, 3, 48) symlat
Definition modmain.f90:344
integer, dimension(maxsymcrys) lsplsymc
Definition modmain.f90:364
Definition modpw.f90:6
integer, dimension(:,:,:), allocatable ihkig
Definition modpw.f90:44
integer, dimension(:,:), allocatable nhk
Definition modpw.f90:40
real(8), dimension(:,:,:,:), allocatable vhkl
Definition modpw.f90:46
pure subroutine r3mtv(a, x, y)
Definition r3mtv.f90:10
subroutine rotaxang(eps, rot, det, v, th)
Definition rotaxang.f90:10