The Elk Code
getevecsv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 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 subroutine getevecsv(fext,ikp,vpl,evecsv)
7 use modmain
8 use modramdisk
9 implicit none
10 ! arguments
11 character(*), intent(in) :: fext
12 integer, intent(in) :: ikp
13 real(8), intent(in) :: vpl(3)
14 complex(8), intent(out) :: evecsv(nstsv,nstsv)
15 ! local variables
16 logical tgs
17 integer isym,lspn
18 integer ik,ist,jst
19 integer recl,nstsv_
20 real(8) vkl_(3),det,v(3),th,t1
21 complex(8) su2(2,2),z1,z2
22 character(256) fname
23 if (ikp > 0) then
24  ik=ikp
25 else
26 ! find the equivalent k-point number and symmetry which rotates vkl to vpl
27  call findkpt(vpl,isym,ik)
28 end if
29 ! construct the filename
30 fname=trim(scrpath)//'EVECSV'//trim(fext)
31 !$OMP CRITICAL(u206)
32 ! read from RAM disk if required
33 if (ramdisk) then
34  call getrd(fname,ik,tgs,v1=vkl_,n1=nstsv_,nzv=nstsv*nstsv,zva=evecsv)
35  if (tgs) goto 10
36 end if
37 ! find the record length
38 inquire(iolength=recl) vkl_,nstsv_,evecsv
39 open(206,file=fname,form='UNFORMATTED',access='DIRECT',recl=recl)
40 read(206,rec=ik) vkl_,nstsv_,evecsv
41 close(206)
42 10 continue
43 !$OMP END CRITICAL(u206)
44 t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
45 if (t1 > epslat) then
46  write(*,*)
47  write(*,'("Error(getevecsv): differing vectors for k-point ",I8)') ik
48  write(*,'(" current : ",3G18.10)') vkl(:,ik)
49  write(*,'(" EVECSV.OUT : ",3G18.10)') vkl_
50  write(*,*)
51  stop
52 end if
53 if (nstsv /= nstsv_) then
54  write(*,*)
55  write(*,'("Error(getevecsv): differing nstsv for k-point ",I8)') ik
56  write(*,'(" current : ",I8)') nstsv
57  write(*,'(" EVECSV.OUT : ",I8)') nstsv_
58  write(*,*)
59  stop
60 end if
61 ! if eigenvectors are spin-unpolarised return
62 if (.not.spinpol) return
63 ! if p = k then return
64 if (ikp > 0) return
65 ! index to global spin rotation in lattice point group
66 lspn=lspnsymc(isym)
67 ! if symmetry element is the identity return
68 if (lspn == 1) return
69 ! find the SU(2) representation of the spin rotation matrix
70 call rotaxang(epslat,symlatc(:,:,lspn),det,v,th)
71 call axangsu2(v,th,su2)
72 ! apply SU(2) matrix to second-variational states (active transformation)
73 do jst=1,nstsv
74  do ist=1,nstfv
75  z1=evecsv(ist,jst)
76  z2=evecsv(ist+nstfv,jst)
77  evecsv(ist,jst)=su2(1,1)*z1+su2(1,2)*z2
78  evecsv(ist+nstfv,jst)=su2(2,1)*z1+su2(2,2)*z2
79  end do
80 end do
81 end subroutine
82 
character(256) scrpath
Definition: modmain.f90:1303
integer, dimension(maxsymcrys) lspnsymc
Definition: modmain.f90:366
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
logical spinpol
Definition: modmain.f90:228
subroutine rotaxang(eps, rot, det, v, th)
Definition: rotaxang.f90:10
logical ramdisk
Definition: modramdisk.f90:9
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
Definition: modramdisk.f90:214
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
pure subroutine axangsu2(v, th, su2)
Definition: axangsu2.f90:10
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24
subroutine findkpt(vpl, isym, ik)
Definition: findkpt.f90:7
integer nstfv
Definition: modmain.f90:887