-
Notifications
You must be signed in to change notification settings - Fork 17
/
ham_bulk.f90
137 lines (105 loc) · 3.3 KB
/
ham_bulk.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
! This subroutine is used to caculate Hamiltonian for
! slab system with open boundary
! History
! Dec/11th/2012 by Quansheng Wu
subroutine ham_bulk(thetax, thetay, hamk_bulk)
use para
implicit none
! wave vector in 2d
real(Dp), intent(in) :: thetax
real(Dp), intent(in) :: thetay
! Hamiltonian of slab system
complex(Dp),intent(out) ::hamk_bulk(ndim, ndim)
! loop index
integer :: i1
integer :: i2
integer :: j1
integer :: j2
complex(Dp) :: ratiox
complex(Dp) :: ratioy
complex(Dp), allocatable :: Hij(:, :, :, :)
allocate(Hij(-ijmax:ijmax, -ijmax:ijmax, nband*2, nband*2))
call ham_qlayer2qlayer3(Hij)
hamk_bulk=0.0d0
do i1= 1, Nx
do j1= 1, Ny
do i2= 1, Nx
do j2= 1, Ny
if (abs(i2-i1).le.ijmax .and. abs(j2-j1).le.ijmax ) then
hamk_bulk( (i1-1)*Nband*2*Ny+(j1-1)*Nband*2+1: &
(i1-1)*Nband*2*Ny+j1*Nband*2, &
(i2-1)*Nband*2*Ny+(j2-1)*Nband*2+1: &
(i2-1)*Nband*2*Ny+j2*Nband*2 ) &
= Hij(i2-i1, j2-j1, 1:nband*2, 1:nband*2)
endif
enddo
enddo
enddo
enddo
return
!>> bulk boundary
ratiox= cos(thetax)+ zi* sin(thetax)
ratioy= cos(thetay)+ zi* sin(thetay)
j1= Ny
j2= 1
do i1= 1, Nx
do i2= 1, Nx
if (abs(i2-i1).le.ijmax)then
hamk_bulk( (i1-1)*Nband*2*Ny+(j1-1)*Nband*2+1: &
(i1-1)*Nband*2*Ny+j1*Nband*2, &
(i2-1)*Nband*2*Ny+(j2-1)*Nband*2+1: &
(i2-1)*Nband*2*Ny+j2*Nband*2 ) &
= Hij(i2-i1, 1, 1:nband*2, 1:nband*2)*ratioy
endif
enddo
enddo
j1= 1
j2= Ny
do i1= 1, Nx
do i2= 1, Nx
if (abs(i2-i1).le.ijmax)then
hamk_bulk( (i1-1)*Nband*2*Ny+(j1-1)*Nband*2+1: &
(i1-1)*Nband*2*Ny+j1*Nband*2, &
(i2-1)*Nband*2*Ny+(j2-1)*Nband*2+1: &
(i2-1)*Nband*2*Ny+j2*Nband*2 ) &
= Hij(i2-i1, -1, 1:nband*2, 1:nband*2)*conjg(ratioy)
endif
enddo
enddo
i1= Nx
i2= 1
do j1= 1, Ny
do j2= 1, Ny
if (abs(j2-j1).le.ijmax)then
hamk_bulk( (i1-1)*Nband*2*Ny+(j1-1)*Nband*2+1: &
(i1-1)*Nband*2*Ny+j1*Nband*2, &
(i2-1)*Nband*2*Ny+(j2-1)*Nband*2+1: &
(i2-1)*Nband*2*Ny+j2*Nband*2 ) &
= Hij(1, j2-j1, 1:nband*2, 1:nband*2)*ratiox
endif
enddo
enddo
i1= 1
i2= Nx
do j1= 1, Ny
do j2= 1, Ny
if (abs(j2-j1).le.ijmax)then
hamk_bulk( (i1-1)*Nband*2*Ny+(j1-1)*Nband*2+1: &
(i1-1)*Nband*2*Ny+j1*Nband*2, &
(i2-1)*Nband*2*Ny+(j2-1)*Nband*2+1: &
(i2-1)*Nband*2*Ny+j2*Nband*2 ) &
= Hij(-1, j2-j1, 1:nband*2, 1:nband*2)*conjg(ratiox)
endif
enddo
enddo
! check hermitcity
do i1=1, Ndim
do i2=1, Ndim
if(abs(hamk_bulk(i1,i2)-conjg(hamk_bulk(i2,i1))).ge.1e-6)then
write(*,*)'there are something wrong with hamk_bulk'
stop
endif
enddo
enddo
return
end