Single-mode approximation (SMA)

1 Introduction

The single-mode approximation (SMA) can be used to analyse excitation spectrum of Quantum dimer models (QDMs).

1.1 Hamilton

H=tr(|||><=|+|=><|||)+vr(|||><|||+|=><=|)

where the first term is the kinetic term, which describes flipping two dimers on every flippable plaquette, and the second term represents a repulsion between two parallel dimers

We focus on a special point, namely, the Rokhsar-Kivelson (RK) point (v=t=1), at which the ground state (|GS>) of the Hamilton is exactly solvable as an equal amplitude superposition of all dimer coverings (|C>).

|GS>=CAC|C>

Dimer covering is that each vertex in the lattice touches exactly one dimer under certain constraints.

1.2 SMA

We define a state |q>=Dα(q)|GS>, where the density operator Dα(q)=1Nrexp(iq.r)nα(r) and nα(r) is dimer number operator. Energy can be defined as:

ω(Q)=12<GS|[D^x(Q),[H,D^x(Q)]]|GS>Dx(Q)

2 code

We present SMA dispersion of the square lattice using fortran language.

2.1 Define parameters

1
2
3
4
5
6
7
8
9
10
11
module configuration
save
integer :: lx=64,nbins=10 !system size:lx*ly
integer :: ly,path
integer :: nn,np !nn:the number of bond along x-direction np: the number of plaquette
integer :: nb !nb: the number of bond alond x and y direction
real(8) :: beta,qx=0,qy=0 !qx,qy: the value of momentum space
integer, allocatable :: bond(:),bondold(:)
integer, allocatable :: bsites(:,:)
complex(8),parameter:: ipi2=(0,1.5707963268),ipi=2*ipi2
end module configuration
1
2
3
4
5
6
7
8
module measurementdata
save
complex :: d_d=0.d0 !<GS|D(-q)D(q)|GS>
complex :: d_hd=0.d0 !<GS|D(-q)HD(q)|GS>
complex :: dh_d=0.d0 !<GS|D(q)HD(-q)|GS>
real(8) :: data1(7)=0.d0
real(8) :: data2(7)=0.d0
end module measurementdata

2.2 Encode each bond of the lattice

The system is made up of np plaquettes, and the bonds of each plaquette are coded in turn

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
subroutine makelattice()
use configuration; implicit none
integer :: s,x,y,z,i
nn=lx*ly
nb=2*nn
np=nn

!|-2-|
!3 4
!|-1-|

allocate(bsites(4,np))
do s=1,nn
bsites(1,s)=s
x=s+lx
if (x>nn) then
bsites(2,s)=x-nn
else
bsites(2,s)=x
endif

bsites(3,s)=s+nn
x=s+nn

if (mod(x,lx)==0) then
bsites(4,s)=x-lx+1
else
bsites(4,s)=x+1
endif
enddo
end subroutine makelattice

2.3 Initial configuration

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
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
!Generate random numbers
real(8) function ran()
implicit none
real(8) :: dmu64
integer(8) :: ran64,mul64,add64
common/bran64/dmu64,ran64,mul64,add64

ran64=ran64*mul64+add64
ran=0.5d0+dmu64*dble(ran64)

end function ran
!---------------------!
subroutine initran(w)
implicit none
integer:: myid, numprocs, namelen, rc,ierr ,comm

integer(8) :: irmax
integer(4) :: w
real(8) :: rdm
real(8) :: dmu64
integer(8) :: ran64,mul64,add64,ii,i
common/bran64/dmu64,ran64,mul64,add64

irmax=2_8**31
irmax=2*(irmax**2-1)+1
mul64=2862933555777941757_8
add64=1013904243
dmu64=0.5d0/dble(irmax)

do ii=1,w
call random_seed()
call random_number(rdm)
enddo
do ii=1,w+3
call random_number(rdm)
enddo
ran64=abs((rdm*mul64)/5+5265361)
end subroutine initran
!-------------------------------------------------------------------!
!Select the initial configuration (dimer covering) based on random numbers
subroutine initconfig()
use configuration; implicit none
integer :: i,b1,b2,b3,b4,ii,next,s1,s2,bs,v0
real(8),external :: ran
real(8) :: rann

allocate(bondold(nb))
allocate(bond(nb))
bond(:)=-1

rann=ran()
print*,rann
if (rann<0.25) then
do i=1,nn
if (mod((i-1)/lx,4)==0) then
if (mod(i-1,4)==0) then
bond(i)=1
bond(i+lx)=1
elseif (mod(i-1,4)==2) then
bond(i+nn)=1
bond(i+nn+1)=1
endif
elseif (mod((i-1)/lx,4)==2) then
if (mod(i-1,4)==2) then
bond(i)=1
bond(i+lx)=1
elseif (mod(i-1,4)==0) then
bond(i+nn)=1
bond(i+nn+1)=1
endif
endif
enddo
elseif (rann<0.5) then
do i=1,nn
if (mod((i-1)/lx,4)==1) then
if (mod(i-1,4)==0) then
bond(i)=1
bond(i+lx)=1
elseif (mod(i-1,4)==2) then
bond(i+nn)=1
bond(i+nn+1)=1
endif
elseif (mod((i-1)/lx,4)==3) then
if (mod(i-1,4)==2) then
bond(i)=1
ii=i+lx
if (ii>nn) ii=ii-nn
bond(ii)=1
elseif (mod(i-1,4)==0) then
bond(i+nn)=1
bond(i+nn+1)=1
endif
endif
enddo
elseif (rann<0.75) then
do i=1,nn
if (mod((i-1)/lx,4)==0) then
if (mod(i-1,4)==1) then
bond(i)=1
bond(i+lx)=1
elseif (mod(i-1,4)==3) then
bond(i+nn)=1
ii=i+nn+1
if (mod(i,lx)==0) ii=ii-lx
bond(ii)=1
endif
elseif (mod((i-1)/lx,4)==2) then
if (mod(i-1,4)==3) then
bond(i)=1
bond(i+lx)=1
elseif (mod(i-1,4)==1) then
bond(i+nn)=1
bond(i+nn+1)=1
endif
endif
enddo
else
do i=1,nn
if (mod((i-1)/lx,4)==1) then
if (mod(i-1,4)==1) then
bond(i)=1
bond(i+lx)=1
elseif (mod(i-1,4)==3) then
bond(i+nn)=1
ii=i+nn+1
if (mod(i,lx)==0) ii=ii-lx
bond(ii)=1
endif
elseif (mod((i-1)/lx,4)==3) then
if (mod(i-1,4)==3) then
bond(i)=1
ii=i+lx
if (ii>nn) ii=ii-nn
bond(ii)=1
elseif (mod(i-1,4)==1) then
bond(i+nn)=1
bond(i+nn+1)=1
endif
endif

enddo
endif
end subroutine initconfig
!-------------------------------------------------------------------!
!any initial configuration can be obtained from another initial configuration by update (|| change to = or = change to ||)
subroutine update()
use configuration; implicit none
integer :: i,ii
real(8),external :: ran
logical,external :: bondcheck
do i=1,np
if (bondcheck(i).and.(ran()<0.5)) then
do ii=1,4
bond(bsites(ii,i))=-bond(bsites(ii,i))
enddo
endif
enddo
end subroutine update

2.4 Calculate

calculate ω ~ k under SMA

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
subroutine measureobservables() 
use configuration; use measurementdata; implicit none
integer :: i,w,x,y,z,ii
complex :: dq,dqh
logical,external :: bondcheck

dq=0
do i=1,nn
x=mod(i-1,lx)
y=(i-1)/lx
dq=dq+exp(ipi*(qx*x+qy*y))*bond(i)
enddo
d_d=d_d+dq*conjg(dq)

do i=1,np
if (bondcheck(i)) then
bondold=bond
do ii=1,4
bondold(bsites(ii,i))=-bondold(bsites(ii,i))
enddo

dqh=0
do ii=1,nn
x=mod(ii-1,lx)
y=(ii-1)/lx
dqh=dqh+exp(ipi*(qx*x+qy*y))*bondold(ii)
enddo

d_hd=-conjg(dqh)*dq+conjg(dq)*dq+d_hd !<GS|D(-q)HD(q)|GS>
dh_d=-dqh*conjg(dq)+dq*conjg(dq)+dh_d !<GS|D(q)HD(-q)|GS>
endif
enddo

end subroutine measureobservables

2.5 main program

we can calculate ω(k) along the Brillouin zone path (0,0)(π,0)(π,π)(0,0)

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
program basic_heisenberg_sse
use configuration;use measurementdata; implicit none

integer :: i,j,msteps=1000,system,ix,iy
character(20)::filename
real(8) :: dqx,dqy

beta=0
ly=lx
dqx=2.d0/lx
dqy=2.d0/ly

call initran(8)
call makelattice()
call initconfig()

path=0

do ix=0,lx/2 ! path 1 from (0,0) to (1,0)
qx=dqx*ix
qy=0

do j=1,nbins
do i=1,msteps
call update()
call measureobservables()
enddo
call writeresults(msteps,j)
enddo
print*,path
path=path+1
enddo

do iy=1,ly/2 ! path 2 from (1,2/ly) to (1,1)
qx=1
qy=dqy*iy
do j=1,nbins
do i=1,msteps
call update()
call measureobservables()
enddo
call writeresults(msteps,j)
enddo
print*,path
path=path+1
enddo

do iy=1,ly/2-1 ! path 3 from (1-2/lx,1-2/ly) to (2/lx,2/ly)
qx=1-dqx*iy
qy=1-dqy*iy
do j=1,nbins
do i=1,msteps
call update()
call measureobservables()
enddo
call writeresults(msteps,j)
enddo
print*,path
path=path+1
enddo

end program basic_heisenberg_sse

2.6 Save

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
subroutine writeresults(msteps,bins)
use configuration; use measurementdata; implicit none

integer :: i,msteps,bins,j
real(8) :: wdata1(7),wdata2(7),t,tt
character(20)::filename

d_d=d_d/msteps
d_hd=d_hd/msteps
dh_d=dh_d/msteps

data1(1)=data1(1)+real(d_hd+dh_d) !molecule
data1(2)=data1(2)+real(d_d) !denominator
data1(3)=data1(3)+real((d_hd+dh_d)/d_d) !$\omega(k)$=molecule/denominator

data2(1)=data2(1)+real(d_hd+dh_d)**2
data2(2)=data2(2)+real(d_d)**2
data2(3)=data2(3)+real((d_hd+dh_d)/d_d)**2

if (bins==nbins) then
do i=1,3
wdata1(i)=data1(i)/bins
wdata2(i)=data2(i)/bins
wdata2(i)=sqrt(abs(wdata2(i)-wdata1(i)**2)/bins) !error
enddo

open(114,file='singlemode.out',position='append')
write(114,*) path,wdata1(3),wdata2(3) !$\omega(k)$,error
close(114)

d_d=0.d0
d_hd=0.d0
dh_d=0.d0
data1=0
data2=0
endif

d_d=0.d0
d_hd=0.d0
dh_d=0.d0

end subroutine writeresults

2.7 result

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
           0   2.3109267291147261E-003   1.6455131249704898E-003
1 6.2533550739288328 2.8371740216404522E-002
2 6.2338413238525394 1.9360616586750323E-002
3 6.1650729179382324 3.1613608837618649E-002
4 6.0500032424926760 1.3714428714990519E-002
5 6.0037383079528812 1.9249298396639518E-002
6 5.9165374279022220 2.1233933146458623E-002
7 5.7899007320404055 2.8688790663723224E-002
8 5.6288352012634277 1.8761953328172567E-002
9 5.4782701969146732 2.0462955105951773E-002
10 5.2886437892913820 1.8608940897921072E-002
11 5.1159195423126222 2.9567157318989505E-002
12 4.9549531936645508 3.1430773070696548E-002
13 4.7688477993011471 2.2452612557261262E-002
14 4.5479550838470457 3.4032958196832191E-002
15 4.3106336593627930 2.7898038113185218E-002
16 4.0874701976776127 1.9144380955909032E-002
17 3.8585180997848512 1.8609396875997147E-002
18 3.6879651308059693 2.8949933454666853E-002
19 3.4377574205398558 1.3578068743458328E-002
20 3.2064547300338746 2.3255134860630623E-002
21 3.0808235883712767 1.7470583451202289E-002
22 2.8069951295852662 2.7647364421329677E-002
23 2.6070729732513427 3.1670710734098995E-002
24 2.3810294151306151 2.9796689027380193E-002
25 2.1672029972076414 1.5971279337898391E-002
26 1.9511149644851684 2.4473548559234566E-002
27 1.7901781082153321 2.7286673718311588E-002
28 1.5865640401840211 3.0836123899655121E-002
29 1.4178456902503966 8.3789057141076002E-003
30 1.2026528716087341 2.2951941343578548E-002
31 1.0390071272850037 3.5212021284243736E-002
32 1.0437846541404725 5.8329403870753194E-002
33 1.1121289670467376 3.9334604709095188E-002
34 1.1652176260948182 3.1558839284690707E-002
35 1.3745294094085694 2.6431431765275463E-002
36 1.4345147013664246 2.1629957788136263E-002
37 1.5875914335250854 2.3702880111387638E-002
38 1.6517681121826171 1.9230563668579375E-002
39 1.7615051031112672 1.4488433511898375E-002
40 1.7648772358894349 1.6782384033197116E-002
41 1.7525395989418029 1.0795956488940417E-002
42 1.7853079319000245 1.5847920523460547E-002
43 1.7468797206878661 7.6477434141953996E-003
44 1.7073505759239196 8.0664708338772963E-003
45 1.6498860597610474 1.0242478225224066E-002
46 1.5684979915618897 8.8734704527437969E-003
47 1.4898424744606018 5.6957569787166405E-003
48 1.3665797352790832 7.4285902427473982E-003
49 1.2612158417701722 8.3717839929171971E-003
50 1.1421930313110351 6.8612852659099147E-003
51 1.0288167774677277 7.4561234512739369E-003
52 0.91099560260772705 5.7387275971854494E-003
53 0.78420571088790891 4.1084180922130774E-003
54 0.66706750988960262 4.9506690982011343E-003
55 0.55218040943145752 5.5888892459151265E-003
56 0.43570946753025053 3.7497197571224258E-003
57 0.34641086459159853 2.3768882259739979E-003
58 0.25858314335346222 3.4487059424092325E-003
59 0.18344641476869583 1.6080926161550495E-003
60 0.11753480210900306 1.8062113973221030E-003
61 6.6618058830499649E-002 5.4123314049076877E-004
62 3.0812994204461576E-002 7.8286098196079113E-004
63 7.0903261657804251E-003 5.8438091546916360E-004
64 0.0000000000000000 0.0000000000000000
65 1.5889596939086915E-002 5.2551528107296812E-004
66 6.1861097440123559E-002 1.4868592619685179E-003
67 0.13374220281839372 1.5077039254814428E-003
!!! do not forget 1/2 (SMA)

References

[1] Läuchli, Andreas M., Sylvain Capponi, and Fakher F. Assaad. “Dynamical dimer correlations at bipartite and non-bipartite Rokhsar–Kivelson points.” Journal of Statistical Mechanics: Theory and Experiment 2008.01 (2008): P01010.

[2] Yan, Zheng, et al. “Widely existing mixed phase structure of the quantum dimer model on a square lattice.” Physical Review B 103.9 (2021): 094421.