The Elk Code
 
Loading...
Searching...
No Matches
dielectric_bse.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 dielectric_bse
7use modmain
8use modomp
9use modtest
10implicit none
11! local variables
12integer a1,a2,ik1,jk1
13integer i1,j1,ist1,jst1
14integer iw,ioc,i,j,nthd
15integer ios,nmbse_
16real(8) e,eji,t1,t2
17complex(8) eta,z1
18character(256) fname
19! allocatable arrays
20real(8), allocatable :: w(:)
21complex(8), allocatable :: pmat(:,:,:),sigma(:,:,:),zv(:,:)
22! initialise global variables
23call init0
24call init1
25! read Fermi energy from a file
26call readfermi
27! get the eigenvalues from file
28do ik1=1,nkpt
29 call getevalsv(filext,ik1,vkl(:,ik1),evalsv(:,ik1))
30end do
31! generate the BSE state index arrays
32call genidxbse
33! allocate global BSE arrays
34if (allocated(evalbse)) deallocate(evalbse)
35allocate(evalbse(nmbse))
36if (allocated(hmlbse)) deallocate(hmlbse)
37allocate(hmlbse(nmbse,nmbse))
38! read in the BSE eigenvectors and eigenvalues
39open(140,file='EVBSE.OUT',form='UNFORMATTED',action='READ',status='OLD', &
40 iostat=ios)
41if (ios /= 0) then
42 write(*,*)
43 write(*,'("Error(dielectric_bse): error opening EVBSE.OUT")')
44 write(*,*)
45 stop
46end if
47read(140) nmbse_
48if (nmbse /= nmbse_) then
49 write(*,*)
50 write(*,'("Error(dielectric_bse): differing nmbse")')
51 write(*,'(" current : ",I6)') nmbse
52 write(*,'(" EVBSE.OUT : ",I6)') nmbse_
53 stop
54end if
55read(140) evalbse
56read(140) hmlbse
57close(140)
58! set up the frequency grid (starting from zero)
59allocate(w(nwplot))
60t1=wplot(2)/dble(nwplot)
61do iw=1,nwplot
62 w(iw)=t1*dble(iw-1)
63end do
64! i divided by the complex relaxation time
65eta=cmplx(0.d0,swidth,8)
66allocate(pmat(nstsv,nstsv,3))
67allocate(sigma(3,3,nwplot))
68allocate(zv(3,nmbse))
69sigma(:,:,:)=0.d0
70zv(:,:)=0.d0
71call holdthd(nmbse,nthd)
72!$OMP PARALLEL DEFAULT(SHARED) &
73!$OMP PRIVATE(ik1,jk1,a1,a2,e,i1,j1) &
74!$OMP PRIVATE(ist1,jst1,eji,z1,i,j) &
75!$OMP REDUCTION(+:sigma) &
76!$OMP NUM_THREADS(nthd)
77! loop over non-reduced k-points
78do ik1=1,nkptnr
79! equivalent reduced k-point
80 jk1=ivkik(ivk(1,ik1),ivk(2,ik1),ivk(3,ik1))
81! read the momentum matrix elements from file
82!$OMP SINGLE
83 call getpmat(vkl(:,ik1),pmat)
84!$OMP END SINGLE
85!$OMP DO SCHEDULE(DYNAMIC)
86 do a2=1,nmbse
87 e=evalbse(a2)
88 do i1=1,nvbse
89 ist1=istbse(i1,ik1)
90 do j1=1,ncbse
91 jst1=jstbse(j1,ik1)
92 a1=ijkbse(i1,j1,ik1)
93 eji=evalsv(jst1,jk1)-evalsv(ist1,jk1)
94 z1=(e/eji)*hmlbse(a1,a2)
95 zv(1:3,a2)=zv(1:3,a2)+z1*pmat(ist1,jst1,1:3)
96 end do
97 end do
98 end do
99!$OMP END DO
100end do
101!$OMP DO SCHEDULE(DYNAMIC)
102do a2=1,nmbse
103 e=evalbse(a2)
104 if (abs(e) > 1.d-8) then
105 do i=1,3
106 do j=1,3
107 z1=zv(i,a2)*conjg(zv(j,a2))/e
108 sigma(i,j,:)=sigma(i,j,:)+z1/(w(:)-e+eta)+conjg(z1)/(w(:)+e+eta)
109 end do
110 end do
111 end if
112end do
113!$OMP END DO
114!$OMP END PARALLEL
115call freethd(nthd)
117sigma(:,:,:)=z1*sigma(:,:,:)
118! loop over tensor components
119do ioc=1,noptcomp
120 i=optcomp(1,ioc)
121 j=optcomp(2,ioc)
122 t1=0.d0
123 if (i == j) t1=1.d0
124 write(fname,'("EPSILON_BSE_",2I1,".OUT")') i,j
125 open(50,file=trim(fname),form='FORMATTED')
126 do iw=1,nwplot
127 t2=t1-fourpi*aimag(sigma(i,j,iw)/(w(iw)+eta))
128 write(50,'(2G18.10)') w(iw),t2
129 end do
130 write(50,*)
131 do iw=1,nwplot
132 t2=fourpi*dble(sigma(i,j,iw)/(w(iw)+eta))
133 write(50,'(2G18.10)') w(iw),t2
134 end do
135 close(50)
136end do
137write(*,*)
138write(*,'("Info(dielectric_bse):")')
139write(*,'(" dielectric tensor written to EPSILON_BSE_ij.OUT")')
140write(*,'(" for components")')
141do ioc=1,noptcomp
142 write(*,'(" i = ",I1,", j = ",I1)') optcomp(1:2,ioc)
143end do
144! write sigma to test file
145call writetest(187,'BSE optical conductivity',nv=nwplot,tol=1.d-2,zva=sigma)
146deallocate(w,pmat,sigma,zv)
147! deallocate global BSE arrays
148deallocate(evalbse,hmlbse)
149end subroutine
150
subroutine dielectric_bse
subroutine genidxbse
Definition genidxbse.f90:7
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine getpmat(vpl, pmat)
Definition getpmat.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
integer nkptnr
Definition modmain.f90:463
integer, dimension(:,:), allocatable istbse
Definition modmain.f90:1192
integer, dimension(3, 27) optcomp
Definition modmain.f90:1090
character(256) filext
Definition modmain.f90:1300
real(8) omega
Definition modmain.f90:20
complex(8), parameter zi
Definition modmain.f90:1239
integer nwplot
Definition modmain.f90:1070
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
complex(8), dimension(:,:), allocatable hmlbse
Definition modmain.f90:1198
real(8), parameter fourpi
Definition modmain.f90:1231
real(8) swidth
Definition modmain.f90:892
integer nkpt
Definition modmain.f90:461
real(8), dimension(2) wplot
Definition modmain.f90:1076
integer, dimension(:,:,:), allocatable ijkbse
Definition modmain.f90:1196
integer, dimension(:,:), allocatable jstbse
Definition modmain.f90:1194
integer noptcomp
Definition modmain.f90:1088
integer nstsv
Definition modmain.f90:886
real(8), dimension(:), allocatable evalbse
Definition modmain.f90:1200
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) occmax
Definition modmain.f90:898
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
integer nmbse
Definition modmain.f90:1190
real(8) wkptnr
Definition modmain.f90:477
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
integer ncbse
Definition modmain.f90:1176
integer nvbse
Definition modmain.f90:1176
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine readfermi
Definition readfermi.f90:10