The Elk Code
getevecfv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007-2008 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: getevecfv
8 ! !INTERFACE:
9 subroutine getevecfv(fext,ikp,vpl,vgpl,evecfv)
10 ! !USES:
11 use modmain
12 use modramdisk
13 ! !INPUT/OUTPUT PARAMETERS:
14 ! fext : filename extension (in,character(*))
15 ! ikp : p-point vector index (in,integer)
16 ! vpl : p-point vector in lattice coordinates (in,real(3))
17 ! vgpl : G+p-vectors in lattice coordinates (out,real(3,ngkmax,nspnfv))
18 ! evecfv : first-variational eigenvectors (out,complex(nmatmax,nstfv,nspnfv))
19 ! !DESCRIPTION:
20 ! Reads in a first-variational eigenvector from file. If the input $k$-point,
21 ! ${\bf p}$, is not in the reduced set, then the eigenvector of the equivalent
22 ! point is read in and the required rotation/translation operations applied.
23 !
24 ! !REVISION HISTORY:
25 ! Created Feburary 2007 (JKD)
26 ! Fixed transformation error, October 2007 (JKD, Anton Kozhevnikov)
27 ! Fixed l.o. rotation, June 2010 (A. Kozhevnikov)
28 !EOP
29 !BOC
30 implicit none
31 ! arguments
32 character(*), intent(in) :: fext
33 integer, intent(in) :: ikp
34 real(8), intent(in) :: vpl(3),vgpl(3,ngkmax,nspnfv)
35 complex(8), intent(out) :: evecfv(nmatmax,nstfv,nspnfv)
36 ! local variables
37 logical tgs
38 integer isym,lspl,ilspl,i,j
39 integer jspn,ist,ilo,l,lm
40 integer ik,ngk_,igp,igk,ig
41 integer is,ia,ja,ias,jas
42 integer recl,nmatmax_,nstfv_,nspnfv_
43 real(8) vkl_(3),v(3)
44 real(8) si(3,3),t1
45 complex(8) z1
46 character(256) fname
47 ! automatic arrays
48 logical done(ngkmax)
49 ! allocatable arrays
50 complex(8), allocatable :: evecfv_(:,:)
51 if (ikp > 0) then
52  ik=ikp
53 else
54 ! find the equivalent k-point number and crystal symmetry element
55  call findkpt(vpl,isym,ik)
56 end if
57 ! construct the filename
58 fname=trim(scrpath)//'EVECFV'//trim(fext)
59 !$OMP CRITICAL(u202)
60 ! read from RAM disk if required
61 if (ramdisk) then
62  call getrd(fname,ik,tgs,v1=vkl_,n1=nmatmax_,n2=nstfv_,n3=nspnfv_, &
63  nzv=nmatmax*nstfv*nspnfv,zva=evecfv)
64  if (tgs) goto 10
65 end if
66 ! find the record length
67 inquire(iolength=recl) vkl_,nmatmax_,nstfv_,nspnfv_,evecfv
68 open(202,file=fname,form='UNFORMATTED',access='DIRECT',recl=recl)
69 read(202,rec=ik) vkl_,nmatmax_,nstfv_,nspnfv_,evecfv
70 close(202)
71 10 continue
72 !$OMP END CRITICAL(u202)
73 t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
74 if (t1 > epslat) then
75  write(*,*)
76  write(*,'("Error(getevecfv): differing vectors for k-point ",I8)') ik
77  write(*,'(" current : ",3G18.10)') vkl(:,ik)
78  write(*,'(" EVECFV.OUT : ",3G18.10)') vkl_
79  write(*,*)
80  stop
81 end if
82 if (nmatmax /= nmatmax_) then
83  write(*,*)
84  write(*,'("Error(getevecfv): differing nmatmax for k-point ",I8)') ik
85  write(*,'(" current : ",I8)') nmatmax
86  write(*,'(" EVECFV.OUT : ",I8)') nmatmax_
87  write(*,*)
88  stop
89 end if
90 if (nstfv /= nstfv_) then
91  write(*,*)
92  write(*,'("Error(getevecfv): differing nstfv for k-point ",I8)') ik
93  write(*,'(" current : ",I8)') nstfv
94  write(*,'(" EVECFV.OUT : ",I8)') nstfv_
95  write(*,*)
96  stop
97 end if
98 if (nspnfv /= nspnfv_) then
99  write(*,*)
100  write(*,'("Error(getevecfv): differing nspnfv for k-point ",I8)') ik
101  write(*,'(" current : ",I8)') nspnfv
102  write(*,'(" EVECFV.OUT : ",I8)') nspnfv_
103  write(*,*)
104  stop
105 end if
106 ! if p = k then return
107 if (ikp > 0) return
108 t1=abs(vpl(1)-vkl(1,ik))+abs(vpl(2)-vkl(2,ik))+abs(vpl(3)-vkl(3,ik))
109 if (t1 < epslat) return
110 ! allocate temporary eigenvector array
111 allocate(evecfv_(nmatmax,nstfv))
112 ! index to spatial rotation in lattice point group
113 lspl=lsplsymc(isym)
114 ! the inverse of the spatial symmetry rotates k into p
115 ilspl=isymlat(lspl)
116 si(1:3,1:3)=dble(symlat(1:3,1:3,ilspl))
117 !-----------------------------------------------!
118 ! translate and rotate APW coefficients !
119 !-----------------------------------------------!
120 ! loop over the first-variational spins
121 do jspn=1,nspnfv
122  ngk_=ngk(jspn,ik)
123  if (tv0symc(isym)) then
124 ! translation vector is zero
125  do ist=1,nstfv
126  do igk=1,ngk_
127  evecfv_(igk,ist)=evecfv(igk,ist,jspn)
128  end do
129  end do
130  else
131 ! non-zero translation vector gives a phase factor
132  v(1:3)=vtcsymc(1:3,isym)
133  do igk=1,ngk_
134  ig=igkig(igk,jspn,ik)
135  t1=vgc(1,ig)*v(1)+vgc(2,ig)*v(2)+vgc(3,ig)*v(3)
136  z1=cmplx(cos(t1),-sin(t1),8)
137  evecfv_(igk,1:nstfv)=z1*evecfv(igk,1:nstfv,jspn)
138  end do
139  end if
140 ! apply spatial rotation operation (passive transformation)
141  done(1:ngk_)=.false.
142  i=1
143  do igk=1,ngk_
144  call r3mtv(si,vgkl(:,igk,jspn,ik),v)
145  j=0
146  do igp=i,ngk_
147  if (done(igp)) cycle
148  if (j == 0) j=igp
149  t1=abs(v(1)-vgpl(1,igp,jspn)) &
150  +abs(v(2)-vgpl(2,igp,jspn)) &
151  +abs(v(3)-vgpl(3,igp,jspn))
152  if (t1 < epslat) then
153  evecfv(igp,1:nstfv,jspn)=evecfv_(igk,1:nstfv)
154  done(igp)=.true.
155  if (igp == j) j=j+1
156  exit
157  end if
158  end do
159  if (j > 0) i=j
160  end do
161 end do
162 !---------------------------------------------------------!
163 ! translate and rotate local-orbital coefficients !
164 !---------------------------------------------------------!
165 if (nlotot > 0) then
166 ! rotate k-point by inverse symmetry matrix
167  call r3mtv(si,vkl(:,ik),v)
168 ! loop over the first-variational spins
169  do jspn=1,nspnfv
170  ngk_=ngk(jspn,ik)
171 ! make a copy of the local-orbital coefficients
172  do i=ngk_+1,nmat(jspn,ik)
173  evecfv_(i,1:nstfv)=evecfv(i,1:nstfv,jspn)
174  end do
175  do is=1,nspecies
176  do ia=1,natoms(is)
177  ias=idxas(ia,is)
178 ! equivalent atom for this symmetry
179  ja=ieqatom(ia,is,isym)
180  jas=idxas(ja,is)
181 ! phase factor from translation
182  t1=-twopi*dot_product(vkl(1:3,ik),atposl(1:3,ja,is))
183  z1=cmplx(cos(t1),sin(t1),8)
184  t1=twopi*dot_product(v(1:3),atposl(1:3,ia,is))
185  z1=z1*cmplx(cos(t1),sin(t1),8)
186 ! rotate local-orbitals (active transformation)
187  do ilo=1,nlorb(is)
188  l=lorbl(ilo,is)
189  lm=l**2+1
190  i=ngk_+idxlo(lm,ilo,ias)
191  j=ngk_+idxlo(lm,ilo,jas)
192  call rotzflm(symlatc(:,:,lspl),l,l,lolmmax,nstfv,nmatmax, &
193  evecfv_(j,1),evecfv(i,1,jspn))
194  evecfv(i:i+2*l,1:nstfv,jspn)=z1*evecfv(i:i+2*l,1:nstfv,jspn)
195  end do
196  end do
197  end do
198  end do
199 end if
200 deallocate(evecfv_)
201 end subroutine
202 !EOC
203 
character(256) scrpath
Definition: modmain.f90:1303
pure subroutine r3mtv(a, x, y)
Definition: r3mtv.f90:10
real(8), parameter twopi
Definition: modmain.f90:1233
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
integer, dimension(:,:,:), allocatable idxlo
Definition: modmain.f90:855
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
logical ramdisk
Definition: modramdisk.f90:9
integer nlotot
Definition: modmain.f90:790
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
integer, dimension(48) isymlat
Definition: modmain.f90:348
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
integer, dimension(:,:,:), allocatable ieqatom
Definition: modmain.f90:368
logical, dimension(maxsymcrys) tv0symc
Definition: modmain.f90:362
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
Definition: modramdisk.f90:214
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
integer, dimension(:,:), allocatable nmat
Definition: modmain.f90:851
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
real(8) epslat
Definition: modmain.f90:24
integer lolmmax
Definition: modmain.f90:800
integer nspecies
Definition: modmain.f90:34
subroutine findkpt(vpl, isym, ik)
Definition: findkpt.f90:7
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
subroutine rotzflm(rot, lmin, lmax, lmmax, n, ld, zflm1, zflm2)
Definition: rotzflm.f90:10
real(8), dimension(3, maxsymcrys) vtcsymc
Definition: modmain.f90:360