The Elk Code
 
Loading...
Searching...
No Matches
getpmat.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 S. Sharma, J. K. Dewhurst 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 getpmat(vpl,pmat)
7use modmain
9implicit none
10! arguments
11real(8), intent(in) :: vpl(3)
12complex(8), intent(out) :: pmat(nstsv,nstsv,3)
13! local variables
14logical tgs
15integer isym,ik,ist,jst
16integer recl,nstsv_
17real(8) vkl_(3),sc(3,3)
18real(8) ei,ej,eij,t1
19complex(8) v1(3),v2(3)
20! find the equivalent k-point number and symmetry which rotates vkl to vpl
21call findkpt(vpl,isym,ik)
22!$OMP CRITICAL(u230)
23! read from RAM disk if required
24if (ramdisk) then
25 call getrd('PMAT.OUT',ik,tgs,v1=vkl_,n1=nstsv_,nzv=nstsv*nstsv*3,zva=pmat)
26 if (tgs) goto 10
27end if
28! find the record length
29inquire(iolength=recl) vkl_,nstsv_,pmat
30open(230,file='PMAT.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
31read(230,rec=ik) vkl_,nstsv_,pmat
32close(230)
3310 continue
34!$OMP END CRITICAL(u230)
35t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
36if (t1 > epslat) then
37 write(*,*)
38 write(*,'("Error(getpmat): differing vectors for k-point ",I8)') ik
39 write(*,'(" current : ",3G18.10)') vkl(:,ik)
40 write(*,'(" PMAT.OUT : ",3G18.10)') vkl_
41 write(*,*)
42 stop
43end if
44if (nstsv /= nstsv_) then
45 write(*,*)
46 write(*,'("Error(getpmat): differing nstsv for k-point ",I8)') ik
47 write(*,'(" current : ",I8)') nstsv
48 write(*,'(" PMAT.OUT : ",I8)') nstsv_
49 write(*,*)
50 stop
51end if
52! apply scissor correction if required
53if (tscissor) then
54 do jst=1,nstsv
55 ej=evalsv(jst,ik)
56 do ist=1,nstsv
57 ei=evalsv(ist,ik)
58 eij=ei-ej
59! note that the eigenvalues have *already* been scissor corrected
60 if ((ei > efermi).and.(ej <= efermi)) then
61 t1=eij/(eij-scissor)
62 else if ((ei <= efermi).and.(ej > efermi)) then
63 t1=eij/(eij+scissor)
64 else
65 cycle
66 end if
67 pmat(ist,jst,1:3)=pmat(ist,jst,1:3)*t1
68 end do
69 end do
70end if
71! if p = k then return
72t1=abs(vpl(1)-vkl(1,ik))+abs(vpl(2)-vkl(2,ik))+abs(vpl(3)-vkl(3,ik))
73if (t1 < epslat) return
74! rotate the matrix elements from the reduced to non-reduced k-point
75sc(1:3,1:3)=symlatc(1:3,1:3,lsplsymc(isym))
76do jst=1,nstsv
77 do ist=1,nstsv
78 v1(1:3)=pmat(ist,jst,1:3)
79 call rz3mv(sc,v1,v2)
80 pmat(ist,jst,1:3)=v2(1:3)
81 end do
82end do
83return
84
85contains
86
87pure subroutine rz3mv(a,x,y)
88implicit none
89real(8), intent(in) :: a(3,3)
90complex(8), intent(in) :: x(3)
91complex(8), intent(out) :: y(3)
92y(1)=a(1,1)*x(1)+a(1,2)*x(2)+a(1,3)*x(3)
93y(2)=a(2,1)*x(1)+a(2,2)*x(2)+a(2,3)*x(3)
94y(3)=a(3,1)*x(1)+a(3,2)*x(2)+a(3,3)*x(3)
95end subroutine
96
97end subroutine
98
subroutine findkpt(vpl, isym, ik)
Definition findkpt.f90:7
subroutine getpmat(vpl, pmat)
Definition getpmat.f90:7
pure subroutine rz3mv(a, x, y)
Definition getpmat.f90:88
logical tscissor
Definition modmain.f90:906
real(8) efermi
Definition modmain.f90:904
real(8), dimension(3, 3, 48) symlatc
Definition modmain.f90:350
real(8) epslat
Definition modmain.f90:24
real(8) scissor
Definition modmain.f90:908
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer, dimension(maxsymcrys) lsplsymc
Definition modmain.f90:364
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
type(file_t), dimension(:), allocatable, private file
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
logical ramdisk
Definition modramdisk.f90:9