1
+ """
2
+ A Tutorial for the Optical SSH Model Ground State Proper
3
+ """
4
+ import numpy
5
+ import h5py
6
+ from renormalizer .model .model import Model , construct_j_matrix
7
+ from renormalizer .model .basis import BasisSimpleElectron , BasisSHO
8
+ from renormalizer .model .op import Op
9
+ from renormalizer .utils import Quantity
10
+ from renormalizer .mps import Mps , Mpo
11
+ from renormalizer .mps .gs import optimize_mps
12
+
13
+
14
+ class OpticalSSHModelGroundState ():
15
+ r"""
16
+ Ground state properties for the optical SSH model.
17
+ Hamiltonian:
18
+ He = sum_{i} t a^\dagger_i a_{i+1}+a_{i+1}^\dagger a_i
19
+ Hph = sum_i w_0 b^\dagger_i b_i
20
+ Heph = sum_{i} g (a^\dagger_i+1 a_i + a^\dagger_i a_{i+1}) (X_{i+1} - X_i)
21
+ X_i = b^\dagger_i + b_i
22
+ """
23
+ def __init__ (self , params ):
24
+ # Model parameters
25
+ self .mol_num = params ['nsites' ]
26
+ self .g = params ['g' ]
27
+ self .w0 = params ['w0' ]
28
+ self .nboson_max = params ['nboson_max' ]
29
+ self .bond_dim = params ['bond_dim' ]
30
+ self .nsweeps = params ['nsweeps' ]
31
+ self .periodic = params ['periodic' ]
32
+ self .t = params ['t' ]
33
+
34
+ # Construct hopping matrix
35
+ j_matrix = construct_j_matrix (self .mol_num , Quantity (self .t ), self .periodic )
36
+
37
+ # Initialize model with basis and Hamiltonian
38
+ self .model = self ._construct_model (j_matrix )
39
+
40
+ def _construct_model (self , j_matrix ):
41
+ # Construct basis
42
+ basis = []
43
+ for imol in range (self .mol_num ):
44
+ basis .append (BasisSimpleElectron (imol ))
45
+ basis .append (BasisSHO ((imol , 0 ), self .w0 , self .nboson_max ))
46
+
47
+ # Construct Hamiltonian terms
48
+ ham = []
49
+ # Electronic hopping terms
50
+ for imol in range (self .mol_num ):
51
+ for jmol in range (self .mol_num ):
52
+ if j_matrix [imol , jmol ] != 0 :
53
+ ham .append (Op (r"a^\dagger a" , [imol , jmol ], j_matrix [imol , jmol ]))
54
+
55
+ # Bosonic terms
56
+ for imol in range (self .mol_num ):
57
+ ham .append (Op ("b^\dagger b" , (imol , 0 ), self .w0 ))
58
+
59
+ # Electron-phonon coupling
60
+ ham .extend (self ._construct_eph_terms ())
61
+
62
+ return Model (basis , ham )
63
+
64
+ def _construct_eph_terms (self ):
65
+ eph_terms = []
66
+ # Bulk terms
67
+ for imol in range (self .mol_num - 1 ):
68
+ eph_terms .extend ([
69
+ Op (r"a^\dagger a" , [imol , imol + 1 ], self .g ) * Op ("b^\dagger+b" , (imol + 1 , 0 )),
70
+ Op (r"a^\dagger a" , [imol , imol + 1 ], - self .g ) * Op ("b^\dagger+b" , (imol , 0 )),
71
+ Op (r"a^\dagger a" , [imol + 1 , imol ], self .g ) * Op ("b^\dagger+b" , (imol + 1 , 0 )),
72
+ Op (r"a^\dagger a" , [imol + 1 , imol ], - self .g ) * Op ("b^\dagger+b" , (imol , 0 ))
73
+ ])
74
+
75
+ # Boundary terms (if periodic)
76
+ if self .periodic :
77
+ eph_terms .extend ([
78
+ Op (r"a^\dagger a" , [self .mol_num - 1 , 0 ], self .g ) * Op ("b^\dagger+b" , (0 , 0 )),
79
+ Op (r"a^\dagger a" , [self .mol_num - 1 , 0 ], - self .g ) * Op ("b^\dagger+b" , (self .mol_num - 1 , 0 )),
80
+ Op (r"a^\dagger a" , [0 , self .mol_num - 1 ], self .g ) * Op ("b^\dagger+b" , (0 , 0 )),
81
+ Op (r"a^\dagger a" , [0 , self .mol_num - 1 ], - self .g ) * Op ("b^\dagger+b" , (self .mol_num - 1 , 0 ))
82
+ ])
83
+ return eph_terms
84
+
85
+ def get_gs_energy (self ):
86
+ # Initialize random MPS and construct MPO
87
+ mps = Mps .random (self .model , 1 , self .bond_dim , percent = 1.0 )
88
+ mpo = Mpo (self .model )
89
+
90
+ # Setup optimization procedure
91
+ procedure = [
92
+ [self .bond_dim // 4 , 0.4 ],
93
+ [self .bond_dim // 2 , 0.2 ],
94
+ [3 * self .bond_dim // 4 , 0.1 ]
95
+ ] + [[self .bond_dim , 0 ]] * (self .nsweeps - 3 )
96
+ mps .optimize_config .procedure = procedure
97
+ mps .optimize_config .method = "2site"
98
+
99
+ # Optimize ground state
100
+ energies , mps = optimize_mps (mps .copy (), mpo )
101
+
102
+ # Calculate observables
103
+ results = {
104
+ 'energies' : energies ,
105
+ 'edof_rdm' : mps .calc_edof_rdm (),
106
+ 'phonon_occupations' : mps .ph_occupations ,
107
+ 'phonon_displacement' : self .calc_phonon_displacement (mps ),
108
+ 'ni_nj' : self .calc_ni_nj (mps )
109
+ }
110
+
111
+ return results
112
+
113
+ def calc_ni_nj (self , mps ):
114
+ """Calculate density-density correlation function."""
115
+ ni_nj = numpy .zeros ((self .mol_num , self .mol_num ))
116
+ for imol in range (self .mol_num ):
117
+ for jmol in range (self .mol_num ):
118
+ ni = Mpo (self .model , Op ("a^\dagger a" , [imol , imol ]))
119
+ nj = Mpo (self .model , Op ("a^\dagger a" , [jmol , jmol ]))
120
+ mpo = ni @ nj
121
+ ni_nj [imol , jmol ] = mps .expectation (mpo )
122
+ return ni_nj
123
+
124
+ def calc_phonon_displacement (self , mps ):
125
+ """Calculate phonon displacement for each site."""
126
+ phonon_displacement = numpy .zeros (self .mol_num )
127
+ for imol in range (self .mol_num ):
128
+ mpo = Mpo (self .model , Op ("b^\dagger+b" , (imol , 0 )))
129
+ phonon_displacement [imol ] = mps .expectation (mpo )
130
+ return phonon_displacement
131
+
132
+
133
+ if __name__ == '__main__' :
134
+ import sys
135
+ # Model parameters
136
+ params = {
137
+ "nsites" : 32 ,
138
+ 'g' : 0.7 ,
139
+ 'w0' : 0.5 ,
140
+ 't' : - 1.0 ,
141
+ 'nboson_max' : 40 ,
142
+ 'bond_dim' : 128 ,
143
+ 'nsweeps' : 20 ,
144
+ 'periodic' : True
145
+ }
146
+
147
+ # Ground state calculation
148
+ job = OpticalSSHModelGroundState (params )
149
+ results = job .get_gs_energy ()
150
+
151
+ # Save results
152
+ with h5py .File ('gs.h5' , 'w' ) as f :
153
+ for key , value in results .items ():
154
+ f .create_dataset (key , data = value )
155
+ f .create_dataset ('gs_energy' , data = min (results ['energies' ]))
0 commit comments