VMPS++
Loading...
Searching...
No Matches
PeierlsHubbardU1xU1.h
Go to the documentation of this file.
1#ifndef HUBBARDMODELU1XU1_H_COMPLEX
2#define HUBBARDMODELU1XU1_H_COMPLEX
3
4#include "symmetry/S1xS2.h"
5#include "symmetry/U1.h"
6#include "bases/FermionBase.h"
8#include "Mpo.h"
9#include "ParamReturner.h"
10#include "Geometry2D.h" // from TOOLS
11
12namespace VMPS
13{
14class PeierlsHubbardU1xU1 : public Mpo<Sym::S1xS2<Sym::U1<Sym::SpinU1>,Sym::U1<Sym::ChargeU1> >,complex<double>>,
15 public HubbardObservables<Sym::S1xS2<Sym::U1<Sym::SpinU1>,Sym::U1<Sym::ChargeU1> >,complex<double>>,
16 public ParamReturner
17{
18public:
19
22 typedef Eigen::Matrix<complex<double>,Eigen::Dynamic,Eigen::Dynamic> MatrixType;
24
25//private:
26
27 typedef Eigen::Index Index;
29
30public:
31
33
34 PeierlsHubbardU1xU1(Mpo<Symmetry,complex<double>> &Mpo_input, const vector<Param> &params)
35 :Mpo<Symmetry,complex<double>>(Mpo_input),
38 {
39 ParamHandler P(params,PeierlsHubbardU1xU1::defaults);
40 size_t Lcell = P.size();
41 N_phys = 0;
42 for (size_t l=0; l<N_sites; ++l) N_phys += P.get<size_t>("Ly",l%Lcell);
43 this->precalc_TwoSiteData();
44 this->HERMITIAN = true;
45 this->HAMILTONIAN = true;
46 };
47
48 PeierlsHubbardU1xU1 (const size_t &L, const vector<Param> &params, const BC &boundary=BC::OPEN, const DMRG::VERBOSITY::OPTION &VERB=DMRG::VERBOSITY::OPTION::ON_EXIT);
49
59 template<typename Symmetry_>
60 static void set_operators (const std::vector<FermionBase<Symmetry_> > &F, const ParamHandler &P,
61 PushType<SiteOperator<Symmetry_,complex<double>>,complex<double>>& pushlist, std::vector<std::vector<std::string>>& labellist,
62 const BC boundary=BC::OPEN);
63
64 static qarray<2> singlet (int N=0) {return qarray<2>{0,N};};
65 static constexpr MODEL_FAMILY FAMILY = HUBBARD;
66 static constexpr int spinfac = 2;
67
68 static const map<string,any> defaults;
69 static const map<string,any> sweep_defaults;
70};
71
72// V is standard next-nearest neighbour density interaction
73// Vz and Vxy are anisotropic isospin-isospin next-nearest neighbour interaction
74const map<string,any> PeierlsHubbardU1xU1::defaults =
75{
76 {"t",1.+0.i}, {"tPrime",0.i}, {"tRung",1.+0.i}, {"tPrimePrime",0.i},
77 {"mu",0.}, {"t0",0.},
78 {"U",0.}, {"Uph",0.},
79 {"V",0.}, {"Vext",0.}, {"Vrung",0.},
80 {"Bz",0.},
81 {"Vz",0.}, {"Vzrung",0.}, {"Vxy",0.}, {"Vxyrung",0.},
82 {"J",0.}, {"Jperp",0.},
83 {"X",0.}, {"Xrung",0.},
84 {"REMOVE_DOUBLE",false}, {"REMOVE_EMPTY",false}, {"REMOVE_UP",false}, {"REMOVE_DN",false}, {"mfactor",1}, {"k",0},
85 {"maxPower",2ul}, {"CYLINDER",false}, {"Ly",1ul}
86};
87
88const map<string,any> PeierlsHubbardU1xU1::sweep_defaults =
89{
90 {"max_alpha",100.}, {"min_alpha",1.}, {"lim_alpha",11ul}, {"eps_svd",1e-7},
91 {"Mincr_abs", 50ul}, {"Mincr_per", 2ul}, {"Mincr_rel", 1.1},
92 {"min_Nsv",0ul}, {"max_Nrich",-1},
93 {"max_halfsweeps",24ul}, {"min_halfsweeps",1ul},
94 {"Minit",2ul}, {"Qinit",2ul}, {"Mlimit",1000ul},
95 {"tol_eigval",1e-7}, {"tol_state",1e-6},
96 {"savePeriod",0ul}, {"CALC_S_ON_EXIT", true}, {"CONVTEST",DMRG::CONVTEST::VAR_2SITE}
97};
98
100PeierlsHubbardU1xU1 (const size_t &L, const vector<Param> &params, const BC &boundary, const DMRG::VERBOSITY::OPTION &VERB)
101:Mpo<Symmetry,complex<double>> (L, Symmetry::qvacuum(), "", PROP::HERMITIAN, PROP::NON_UNITARY, boundary, VERB),
102 HubbardObservables(L,params,PeierlsHubbardU1xU1::defaults),
103 ParamReturner(PeierlsHubbardU1xU1::sweep_defaults)
104{
105 ParamHandler P(params,defaults);
106 size_t Lcell = P.size();
107
108 for (size_t l=0; l<N_sites; ++l)
109 {
110 N_phys += P.get<size_t>("Ly",l%Lcell);
111 setLocBasis(F[l].get_basis().qloc(),l);
112 }
113
114 this->set_name("Peierls-Hubbard");
115
116 PushType<SiteOperator<Symmetry,complex<double>>,complex<double>> pushlist;
117 std::vector<std::vector<std::string>> labellist;
118 set_operators(F, P, pushlist, labellist, boundary);
119
120 this->construct_from_pushlist(pushlist, labellist, Lcell);
121 this->finalize(PROP::COMPRESS, P.get<size_t>("maxPower"));
122
123 this->precalc_TwoSiteData();
124}
125
126template<typename Symmetry_>
128set_operators (const std::vector<FermionBase<Symmetry_> > &F, const ParamHandler &P, PushType<SiteOperator<Symmetry_,complex<double>>,complex<double>>& pushlist, std::vector<std::vector<std::string>>& labellist, const BC boundary)
129{
130 std::size_t Lcell = P.size();
131 std::size_t N_sites = F.size();
132 if(labellist.size() != N_sites) {labellist.resize(N_sites);}
133
134 for (std::size_t loc=0; loc<N_sites; ++loc)
135 {
136 size_t lp1 = (loc+1)%N_sites;
137 size_t lp2 = (loc+2)%N_sites;
138 size_t lp3 = (loc+3)%N_sites;
139
140 std::size_t orbitals = F[loc].orbitals();
141 std::size_t next_orbitals = F[lp1].orbitals();
142 std::size_t nextn_orbitals = F[lp2].orbitals();
143 std::size_t nnextn_orbitals = F[lp3].orbitals();
144
145 stringstream ss;
146 ss << "Ly=" << P.get<size_t>("Ly",loc%Lcell);
147 labellist[loc].push_back(ss.str());
148
149 auto push_full = [&N_sites, &loc, &F, &P, &pushlist, &labellist, &boundary] (string xxxFull, string label,
150 const vector<SiteOperatorQ<Symmetry_,MatrixType>> &first,
151 const vector<vector<SiteOperatorQ<Symmetry_,MatrixType>>> &last,
152 vector<complex<double>> factor,
153 vector<bool> CONJ,
154 bool FERMIONIC) -> void
155 {
156 ArrayXXcd Full = P.get<Eigen::ArrayXXcd>(xxxFull);
157 vector<vector<std::pair<size_t,complex<double>> > > R = Geometry2D::rangeFormat(Full);
158
159 if (static_cast<bool>(boundary)) {assert(R.size() == N_sites and "Use an (N_sites)x(N_sites) hopping matrix for open BC!");}
160 else {assert(R.size() >= 2*N_sites and "Use at least a (2*N_sites)x(N_sites) hopping matrix for infinite BC!");}
161
162 for (size_t j=0; j<first.size(); j++)
163 for (size_t h=0; h<R[loc].size(); ++h)
164 {
165 size_t range = R[loc][h].first;
166 complex<double> value = R[loc][h].second;
167
168 if (range != 0)
169 {
170 vector<SiteOperatorQ<Symmetry_,MatrixType> > ops(range+1);
171 ops[0] = first[j];
172 for (size_t i=1; i<range; ++i)
173 {
174 if (FERMIONIC) {ops[i] = F[(loc+i)%N_sites].sign().template cast<complex<double>>();}
175 else {ops[i] = F[(loc+i)%N_sites].Id().template cast<complex<double>>();}
176 }
177 ops[range] = last[j][(loc+range)%N_sites];
178 complex<double> total_value = factor[j] * value;
179 if (CONJ[j]) total_value = conj(total_value);
180 pushlist.push_back(std::make_tuple(loc, ops, total_value));
181 }
182 }
183
184 stringstream ss;
185 ss << label << "(" << Geometry2D::hoppingInfo(Full) << ")";
186 labellist[loc].push_back(ss.str());
187 };
188
189 if (P.HAS("tFull"))
190 {
191 SiteOperatorQ<Symmetry_,MatrixType> cUP_sign_local = (F[loc].c(UP,0) * F[loc].sign()).template cast<complex<double>>();
192 SiteOperatorQ<Symmetry_,MatrixType> cDN_sign_local = (F[loc].c(DN,0) * F[loc].sign()).template cast<complex<double>>();
193 SiteOperatorQ<Symmetry_,MatrixType> cdagUP_sign_local = (F[loc].cdag(UP,0) * F[loc].sign()).template cast<complex<double>>();
194 SiteOperatorQ<Symmetry_,MatrixType> cdagDN_sign_local = (F[loc].cdag(DN,0) * F[loc].sign()).template cast<complex<double>>();
195
196 vector<SiteOperatorQ<Symmetry_,MatrixType> > cUP_ranges(N_sites); for (size_t i=0; i<N_sites; ++i) {cUP_ranges[i] = F[i].c(UP,0).template cast<complex<double>>();}
197 vector<SiteOperatorQ<Symmetry_,MatrixType> > cdagUP_ranges(N_sites); for (size_t i=0; i<N_sites; ++i) {cdagUP_ranges[i] = F[i].cdag(UP,0).template cast<complex<double>>();}
198 vector<SiteOperatorQ<Symmetry_,MatrixType> > cDN_ranges(N_sites); for (size_t i=0; i<N_sites; ++i) {cDN_ranges[i] = F[i].c(DN,0).template cast<complex<double>>();}
199 vector<SiteOperatorQ<Symmetry_,MatrixType> > cdagDN_ranges(N_sites); for (size_t i=0; i<N_sites; ++i) {cdagDN_ranges[i] = F[i].cdag(DN,0).template cast<complex<double>>();}
200
201 vector<SiteOperatorQ<Symmetry_,MatrixType> > frst {cdagUP_sign_local, cUP_sign_local, cdagDN_sign_local, cDN_sign_local};
202 vector<vector<SiteOperatorQ<Symmetry_,MatrixType> > > last {cUP_ranges, cdagUP_ranges, cDN_ranges, cdagDN_ranges};
203 push_full("tFull", "tᵢⱼ", frst, last, {-1., +1., -1., +1.}, {false, true, false, true}, PROP::FERMIONIC);
204 }
205
206 if (P.HAS("Vzfull"))
207 {
208 vector<SiteOperatorQ<Symmetry_,MatrixType> > first {F[loc].Tz(0).template cast<complex<double>>()};
209 vector<SiteOperatorQ<Symmetry_,MatrixType> > Tz_ranges(N_sites);
210 for (size_t i=0; i<N_sites; i++)
211 {
212 Tz_ranges[i] = F[i].Tz(0).template cast<complex<double>>();
213 }
214
215 vector<vector<SiteOperatorQ<Symmetry_,MatrixType> > > last {Tz_ranges};
216 push_full("Vzfull", "Vzᵢⱼ", first, last, {1.}, {false}, PROP::BOSONIC);
217 }
218
219 if (P.HAS("vzfull"))
220 {
221 vector<SiteOperatorQ<Symmetry_,MatrixType> > first {F[loc].tz(0).template cast<complex<double>>()};
222 vector<SiteOperatorQ<Symmetry_,MatrixType> > tz_ranges(N_sites);
223 for (size_t i=0; i<N_sites; i++)
224 {
225 tz_ranges[i] = F[i].tz(0).template cast<complex<double>>();
226 }
227
228 vector<vector<SiteOperatorQ<Symmetry_,MatrixType> > > last {tz_ranges};
229 push_full("vzfull", "vzᵢⱼ", first, last, {1.}, {false}, PROP::BOSONIC);
230 }
231
232 if (P.HAS("VextFull"))
233 {
234 vector<SiteOperatorQ<Symmetry_,MatrixType> > first {F[loc].n(0).template cast<complex<double>>()};
235 vector<SiteOperatorQ<Symmetry_,MatrixType> > n_ranges(N_sites); for (size_t i=0; i<N_sites; i++) {n_ranges[i] = F[i].n(0).template cast<complex<double>>();}
236 vector<vector<SiteOperatorQ<Symmetry_,MatrixType> > > last {n_ranges};
237 push_full("VextFull", "Vextᵢⱼ", first, last, {1.}, {false}, PROP::BOSONIC);
238 }
239
240 if (P.HAS("Jfull"))
241 {
242 vector<SiteOperatorQ<Symmetry_,MatrixType> > first {F[loc].Sp(0).template cast<complex<double>>(),
243 F[loc].Sm(0).template cast<complex<double>>(),
244 F[loc].Sz(0).template cast<complex<double>>()};
245 vector<SiteOperatorQ<Symmetry_,MatrixType> > Sp_ranges(N_sites);
246 vector<SiteOperatorQ<Symmetry_,MatrixType> > Sm_ranges(N_sites);
247 vector<SiteOperatorQ<Symmetry_,MatrixType> > Sz_ranges(N_sites);
248 for (size_t i=0; i<N_sites; i++)
249 {
250 Sp_ranges[i] = F[i].Sp(0).template cast<complex<double>>();
251 Sm_ranges[i] = F[i].Sm(0).template cast<complex<double>>();
252 Sz_ranges[i] = F[i].Sz(0).template cast<complex<double>>();
253 }
254
255 vector<vector<SiteOperatorQ<Symmetry_,MatrixType> > > last {Sm_ranges, Sp_ranges, Sz_ranges};
256 push_full("Jfull", "Jᵢⱼ", first, last, {0.5,0.5,1.}, {false,false,false}, PROP::BOSONIC);
257 }
258
259 // Local terms: U, t0, μ, t⟂, V⟂, J⟂
260 param1d U = P.fill_array1d<double>("U", "Uorb", orbitals, loc%Lcell);
261 param1d Uph = P.fill_array1d<double>("Uph", "Uphorb", orbitals, loc%Lcell);
262 param1d t0 = P.fill_array1d<double>("t0", "t0orb", orbitals, loc%Lcell);
263 param1d mu = P.fill_array1d<double>("mu", "muorb", orbitals, loc%Lcell);
264 param1d Bz = P.fill_array1d<double>("Bz", "Bzorb", orbitals, loc%Lcell);
265 param2d tperp = P.fill_array2d<complex<double>>("tRung", "t", "tPerp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
266 param2d Vperp = P.fill_array2d<double>("VRung", "V", "VPerp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
267 param2d Vzperp = P.fill_array2d<double>("VzRung", "Vz", "VzPerp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
268 param2d Vxyperp = P.fill_array2d<double>("VxyRung", "Vxy", "VxyPerp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
269 param2d Jperp = P.fill_array2d<double>("JRung", "J", "JPerp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
270
271 labellist[loc].push_back(U.label);
272 labellist[loc].push_back(Uph.label);
273 labellist[loc].push_back(t0.label);
274 labellist[loc].push_back(mu.label);
275 labellist[loc].push_back(Bz.label);
276 labellist[loc].push_back(tperp.label);
277 labellist[loc].push_back(Vperp.label);
278 labellist[loc].push_back(Vzperp.label);
279 labellist[loc].push_back(Vxyperp.label);
280 labellist[loc].push_back(Jperp.label);
281 ArrayXd C_array = F[loc].ZeroField();
282
284 (
285 F[loc].template HubbardHamiltonian<complex<double>,Symmetry_>(U.a.cast<complex<double>>(),
286 Uph.a.cast<complex<double>>(),
287 (t0.a-mu.a).cast<complex<double>>(),
288 Bz.a.cast<complex<double>>(),
289 tperp.a,
290 Vperp.a.cast<complex<double>>(),
291 Vzperp.a.cast<complex<double>>(),
292 Vxyperp.a.cast<complex<double>>(),
293 Jperp.a.cast<complex<double>>(),
294 Jperp.a.cast<complex<double>>(),
295 C_array.cast<complex<double>>())
296 );
297 pushlist.push_back(std::make_tuple(loc, Hloc, 1.+0.i));
298 }
299}
300
301} // end namespace VMPS::models
302
303#endif
@ DN
Definition: DmrgTypedefs.h:38
@ UP
Definition: DmrgTypedefs.h:37
MODEL_FAMILY
Definition: DmrgTypedefs.h:96
@ HUBBARD
Definition: DmrgTypedefs.h:96
BC
Definition: DmrgTypedefs.h:161
std::enable_if< Dummy::IS_SPIN_SU2() and!Dummy::IS_CHARGE_SU2(), Mpo< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::U1< Sym::ChargeU1 > >, complex< double > > >::type P(size_t locx1, size_t locx2, size_t locy1=0, size_t locy2=0) const
std::enable_if<!Dummy::IS_SPIN_SU2(), Mpo< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::U1< Sym::ChargeU1 > >, complex< double > > >::type Sz(size_t locx, size_t locy=0) const
Mpo< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::U1< Sym::ChargeU1 > >, complex< double > > Id() const
std::enable_if<!Dummy::IS_SPIN_SU2(), Mpo< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::U1< Sym::ChargeU1 > >, complex< double > > >::type Sm(size_t locx, size_t locy=0) const
vector< FermionBase< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::U1< Sym::ChargeU1 > > > > F
std::size_t N_phys
Definition: MpoTerms.h:400
void finalize(const bool COMPRESS=true, const std::size_t power=1, const double tolerance=::mynumeric_limits< double >::epsilon())
Definition: MpoTerms.h:1281
std::size_t N_sites
Definition: MpoTerms.h:395
void setLocBasis(const std::vector< std::vector< qType > > &q)
Definition: MpoTerms.h:715
DMRG::VERBOSITY::OPTION VERB
Definition: MpoTerms.h:102
MpoTerms< Symmetry, OtherScalar > cast()
Definition: MpoTerms.h:3076
void set_name(const std::string &label_in)
Definition: MpoTerms.h:471
std::string label
Definition: MpoTerms.h:615
Definition: Mpo.h:40
static std::vector< T > get_N_site_interaction(T const &Op0, Operator const &... Ops)
Definition: Mpo.h:117
void construct_from_pushlist(const PushType< OperatorType, CouplScalar > &pushlist, const std::vector< std::vector< std::string > > &labellist, size_t Lcell)
Definition: U1.h:25
static constexpr int spinfac
static constexpr MODEL_FAMILY FAMILY
PeierlsHubbardU1xU1(Mpo< Symmetry, complex< double > > &Mpo_input, const vector< Param > &params)
static void set_operators(const std::vector< FermionBase< Symmetry_ > > &F, const ParamHandler &P, PushType< SiteOperator< Symmetry_, complex< double > >, complex< double > > &pushlist, std::vector< std::vector< std::string > > &labellist, const BC boundary=BC::OPEN)
Eigen::Matrix< complex< double >, Eigen::Dynamic, Eigen::Dynamic > MatrixType
static qarray< 2 > singlet(int N=0)
static const map< string, any > defaults
static const map< string, any > sweep_defaults
Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::U1< Sym::ChargeU1 > > Symmetry
#define MAKE_TYPEDEFS(MODEL)
Definition: macros.h:4
const bool COMPRESS
Definition: DmrgTypedefs.h:499
const bool FERMIONIC
Definition: DmrgTypedefs.h:496
const bool BOSONIC
Definition: DmrgTypedefs.h:498
Definition: qarray.h:26