The Elk Code
 
Loading...
Searching...
No Matches
hmlephb.f90
Go to the documentation of this file.
1
2! Copyright (C) 2020 Chung-Yu Wang, J. K. Dewhurst, S. Sharma and
3! E. K. U. Gross. This file is distributed under the terms of the GNU General
4! Public License. See the file COPYING for license details.
5
6subroutine hmlephb(ik,b)
7use modmain
8use modphonon
9use modbog
10use modomp
11implicit none
12! arguments
13integer, intent(in) :: ik
14complex(8), intent(out) :: b(nstsv,nstsv)
15! local variables
16integer iq,jq,ikq,isym,nthd
17integer i1,i2,j1,j2,l1,l2
18real(8) vl(3),t0
19complex(8) z1,z2
20! automatic arrays
21complex(4) ephmat(nstsv,nstsv,nbph)
22complex(8) x(nbph,nstsv),y(nstsv,nbph)
23! prefactor
24t0=-2.d0*wqptnr*ephscf(1)**2/dengy
25b(:,:)=0.d0
26if (.not.anomalous) return
27! parallel loop over non-reduced q-points
28call holdthd(nqptnr,nthd)
29!$OMP PARALLEL DO DEFAULT(SHARED) &
30!$OMP PRIVATE(ephmat,x,y,jq,vl,isym,ikq) &
31!$OMP PRIVATE(i1,i2,j1,j2,l1,l2,z1,z2) &
32!$OMP REDUCTION(+:b) &
33!$OMP SCHEDULE(DYNAMIC) &
34!$OMP NUM_THREADS(nthd)
35do iq=1,nqptnr
36! equivalent reduced q-point
37 jq=ivqiq(ivq(1,iq),ivq(2,iq),ivq(3,iq))
38! k+q-vector in lattice coordinates
39 vl(:)=vkl(:,ik)+vql(:,iq)
40! index to reduced k+q-vector
41 call findkpt(vl,isym,ikq)
42! read in the electron-phonon matrix elements from file
43 call getephmkq(iq,ik,ephmat)
44! perform the contraction
45 do i2=1,nstsv
46 if (abs(evalsv(i2,ik)-efermi) > ecutb) cycle
47 do j1=1,nstsv
48 do l2=1,nbph
49 z1=0.d0
50 do j2=1,nstsv
51 z1=z1+ephmat(j2,i2,l2)*duv(j1,j2,ikq)
52 end do
53 x(l2,j1)=z1
54 end do
55 end do
56 do l1=1,nbph
57 do j1=1,nstsv
58 z1=0.d0
59 do l2=1,nbph
60 z2=dxx(l2,l1,jq)+dwx(l2,l1,jq)
61 z1=z1+z2*x(l2,j1)
62 end do
63 y(j1,l1)=z1
64 end do
65 end do
66 do i1=1,nstsv
67 if (bdiag.and.(i1 /= i2)) cycle
68 if (abs(evalsv(i1,ik)-efermi) > ecutb) cycle
69 z1=0.d0
70 do l1=1,nbph
71 do j1=1,nstsv
72 z1=z1+conjg(ephmat(j1,i1,l1))*y(j1,l1)
73 end do
74 end do
75 b(i1,i2)=b(i1,i2)+t0*z1
76 end do
77 end do
78end do
79!$OMP END PARALLEL DO
80call freethd(nthd)
81end subroutine
82
subroutine findkpt(vpl, isym, ik)
Definition findkpt.f90:7
subroutine getephmkq(iqp, ikp, ephmat)
Definition getephmkq.f90:7
subroutine hmlephb(ik, b)
Definition hmlephb.f90:7
complex(8), dimension(:,:,:), pointer, contiguous duv
Definition modbog.f90:25
logical bdiag
Definition modbog.f90:29
complex(8), dimension(:,:,:), pointer, contiguous dxx
Definition modbog.f90:43
real(8) ecutb
Definition modbog.f90:31
complex(8), dimension(:,:,:), pointer, contiguous dwx
Definition modbog.f90:43
real(8) efermi
Definition modmain.f90:904
integer, dimension(:,:), allocatable ivq
Definition modmain.f90:529
integer nqptnr
Definition modmain.f90:527
integer, dimension(:,:,:), allocatable ivqiq
Definition modmain.f90:531
real(8), dimension(:,:), allocatable vql
Definition modmain.f90:545
real(8) wqptnr
Definition modmain.f90:551
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
logical anomalous
real(8), dimension(2) ephscf
real(8) dengy