VMPS++
Loading...
Searching...
No Matches
HubbardU1xSU2.h
Go to the documentation of this file.
1#ifndef HUBBARDMODEL_U1XSU2_H_
2#define HUBBARDMODEL_U1XSU2_H_
3
4#include "symmetry/S1xS2.h"
5#include "symmetry/U1.h"
6#include "symmetry/SU2.h"
8#include "bases/FermionBase.h"
10#include "Mpo.h"
11#include "ParamReturner.h"
12#include "Geometry2D.h" // from TOOLS
13
14namespace VMPS
15{
16class HubbardU1xSU2 : public Mpo<Sym::S1xS2<Sym::U1<Sym::SpinU1>,Sym::SU2<Sym::ChargeSU2> > ,double>,
17 public HubbardObservables<Sym::S1xS2<Sym::U1<Sym::SpinU1>,Sym::SU2<Sym::ChargeSU2> > >,
18 public ParamReturner
19{
20public:
23
24private:
25
26 typedef Eigen::Index Index;
28 typedef Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic> MatrixType;
29
31
32public:
33
36
37 HubbardU1xSU2(Mpo<Symmetry> &Mpo_input, const vector<Param> &params)
38 :Mpo<Symmetry>(Mpo_input),
41 {
42 ParamHandler P(params,HubbardU1xSU2::defaults);
43 size_t Lcell = P.size();
44 N_phys = 0;
45 for (size_t l=0; l<N_sites; ++l) N_phys += P.get<size_t>("Ly",l%Lcell);
46 this->precalc_TwoSiteData();
47 this->HERMITIAN = true;
48 this->HAMILTONIAN = true;
49 };
50
51 HubbardU1xSU2 (const size_t &L, const vector<Param> &params, const BC &boundary=BC::OPEN, const DMRG::VERBOSITY::OPTION &VERB=DMRG::VERBOSITY::OPTION::ON_EXIT);
53
63 template<typename Symmetry_>
64 static void set_operators (const std::vector<FermionBase<Symmetry_> > &F,
65 const vector<SUB_LATTICE> &G, const ParamHandler &P,
66 PushType<SiteOperator<Symmetry_,double>,double>& pushlist, std::vector<std::vector<std::string>>& labellist,
67 const BC boundary=BC::OPEN);
68
69 static qarray<2> singlet (int N=0, int L=0)
70 {
71 assert(N%2==0);
72 int T = abs(0.5*(N-L));
73 return qarray<2>{0,2*T+1};
74 };
75 static constexpr MODEL_FAMILY FAMILY = HUBBARD;
76 static constexpr int spinfac = 2;
77
78 static const map<string,any> defaults;
79 static const map<string,any> sweep_defaults;
80};
81
82const map<string,any> HubbardU1xSU2::defaults =
83{
84 {"t",1.}, {"tRung",1.}, {"tPrimePrime",0.},
85 {"Uph",0.},
86 {"V",0.}, {"Vrung",0.},
87 {"J",0.}, {"Jrung",0.},
88 {"X",0.}, {"Xrung",0.},
89 {"REMOVE_DOUBLE",false}, {"REMOVE_EMPTY",false}, {"REMOVE_UP",false}, {"REMOVE_DN",false}, {"mfactor",1}, {"k",1},
90 {"maxPower",2ul}, {"CYLINDER",false}, {"Ly",1ul}
91};
92
93const map<string,any> HubbardU1xSU2::sweep_defaults =
94{
95 {"max_alpha",100.}, {"min_alpha",1e-11}, {"lim_alpha",11ul}, {"eps_svd",1e-7},
96 {"Dincr_abs", 2ul}, {"Dincr_per", 2ul}, {"Dincr_rel", 1.1},
97 {"min_Nsv",0ul}, {"max_Nrich",-1},
98 {"max_halfsweeps",30ul}, {"min_halfsweeps",6ul},
99 {"Dinit",4ul}, {"Qinit",10ul}, {"Dlimit",500ul},
100 {"tol_eigval",1e-6}, {"tol_state",1e-5},
101 {"savePeriod",0ul}, {"CALC_S_ON_EXIT", true}, {"CONVTEST", DMRG::CONVTEST::VAR_2SITE}
102};
103
105HubbardU1xSU2 (const size_t &L, const vector<Param> &params, const BC &boundary, const DMRG::VERBOSITY::OPTION &VERB)
106:Mpo<Symmetry> (L, qarray<Symmetry::Nq>({0,1}), "", PROP::HERMITIAN, PROP::NON_UNITARY, boundary, VERB),
109{
110 ParamHandler P(params,defaults);
111 size_t Lcell = P.size();
112
113// assert(Lcell > 1 and "You need to set a unit cell with at least Lcell=2 for the charge-SU(2) symmetry!");
114 for (size_t l=0; l<N_sites; ++l)
115 {
116 N_phys += P.get<size_t>("Ly",l%Lcell);
117 setLocBasis(F[l].get_basis().qloc(),l);
118 }
119
120 this->set_name("Hubbard");
121
123 std::vector<std::vector<std::string>> labellist;
124 set_operators(F, G, P, pushlist, labellist, boundary); // F, G are set in HubbardObservables
125
126 this->construct_from_pushlist(pushlist, labellist, Lcell);
127 this->finalize(PROP::COMPRESS, P.get<size_t>("maxPower"));
128
129 this->precalc_TwoSiteData();
130}
131
132template<typename Symmetry_>
134set_operators (const std::vector<FermionBase<Symmetry_> > &F, const vector<SUB_LATTICE> &G, const ParamHandler &P, PushType<SiteOperator<Symmetry_,double>,double>& pushlist, std::vector<std::vector<std::string>>& labellist, const BC boundary)
135{
136 std::size_t Lcell = P.size();
137 std::size_t N_sites = F.size();
138 if(labellist.size() != N_sites) {labellist.resize(N_sites);}
139
140 for(std::size_t loc=0; loc<N_sites; ++loc)
141 {
142 size_t lp1 = (loc+1)%N_sites;
143 size_t lp2 = (loc+2)%N_sites;
144 size_t lp3 = (loc+3)%N_sites;
145
146 //auto Gloc = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,loc)));
147 //auto Glp1 = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,lp1)));
148 //auto Glp2 = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,lp2)));
149 //auto Glp3 = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,lp3)));
150 //lout << G[loc] << "\t" << G[lp1] << "\t" << G[lp2] << "\t" << G[lp3] << endl;
151
152 std::size_t orbitals = F[loc].orbitals();
153 std::size_t next_orbitals = F[lp1].orbitals();
154 std::size_t next3_orbitals = F[lp3].orbitals();
155
156 stringstream ss;
157 ss << "Ly=" << P.get<size_t>("Ly",loc%Lcell);
158 labellist[loc].push_back(ss.str());
159
160 auto push_full = [&N_sites, &loc, &F, &P, &pushlist, &labellist, &boundary] (string xxxFull, string label,
161 const vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > &first,
162 const vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > &last,
163 vector<double> factor, bool FERMIONIC) -> void
164 {
165 ArrayXXd Full = P.get<Eigen::ArrayXXd>(xxxFull);
166 vector<vector<std::pair<size_t,double> > > R = Geometry2D::rangeFormat(Full);
167
168 if (static_cast<bool>(boundary)) {assert(R.size() == N_sites and "Use an (N_sites)x(N_sites) hopping matrix for open BC!");}
169 else {assert(R.size() >= 2*N_sites and "Use at least a (2*N_sites)x(N_sites) hopping matrix for infinite BC!");}
170
171 for (size_t j=0; j<first.size(); j++)
172 for (size_t h=0; h<R[loc].size(); ++h)
173 {
174 size_t range = R[loc][h].first;
175 double value = R[loc][h].second;
176
177 if (range != 0)
178 {
179 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > ops(range+1);
180 ops[0] = first[j];
181 for (size_t i=1; i<range; ++i)
182 {
183 if (FERMIONIC) {ops[i] = F[(loc+i)%N_sites].sign();}
184 else {ops[i] = F[(loc+i)%N_sites].Id();}
185 }
186 ops[range] = last[j][(loc+range)%N_sites];
187 pushlist.push_back(std::make_tuple(loc, ops, factor[j] * value));
188 }
189 }
190
191 stringstream ss;
192 ss << label << "(" << Geometry2D::hoppingInfo(Full) << ")";
193 labellist[loc].push_back(ss.str());
194 };
195
196 if (P.HAS("tFull"))
197 {
198 SiteOperatorQ<Symmetry_,Eigen::MatrixXd> cdagup_sign_local = F[loc].cdag(UP,G[loc],0) * F[loc].sign();
199 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > cup_ranges(N_sites);
200 SiteOperatorQ<Symmetry_,Eigen::MatrixXd> cdagdn_sign_local = F[loc].cdag(DN,G[loc],0) * F[loc].sign();
201 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > cdn_ranges(N_sites);
202 for (size_t i=0; i<N_sites; i++)
203 {
204 //auto Gi = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,i)));
205 cup_ranges[i] = F[i].c(UP,G[i],0);
206 }
207 for (size_t i=0; i<N_sites; i++)
208 {
209 //auto Gi = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,i)));
210 cdn_ranges[i] = F[i].c(DN,G[i],0);
211 }
212
213 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > frst {cdagup_sign_local, cdagdn_sign_local};
214 vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > last {cup_ranges, cdn_ranges};
215 push_full("tFull", "tᵢⱼ", frst, last, {-std::sqrt(2.),-std::sqrt(2.)}, PROP::FERMIONIC);
216 }
217
218 // not implemented for this symmetry:
219// if (P.HAS("Vfull"))
220// {
221// vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > first {F[loc].Tdag(0)};
222// vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > T_ranges(N_sites); for (size_t i=0; i<N_sites; i++) {T_ranges[i] = F[i].T(0);}
223// vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > last {T_ranges};
224// push_full("Vfull", "Jᵢⱼ", first, last, {std::sqrt(3.)}, PROP::BOSONIC);
225// }
226//
227// if (P.HAS("Jfull"))
228// {
229// vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > first {F[loc].Sdag(0)};
230// vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > S_ranges(N_sites); for (size_t i=0; i<N_sites; i++) {S_ranges[i] = F[i].S(0);}
231// vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > last {S_ranges};
232// push_full("Jfull", "Jᵢⱼ", first, last, {std::sqrt(3.)}, PROP::BOSONIC);
233// }
234//
235// if (P.HAS("Xfull"))
236// {
237// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsidagRloc = ((F[loc].ns() * F[loc].cdag(G[loc])) * F[loc].sign());
238// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsidagLloc = ((F[loc].cdag(G[loc]) * F[loc].sign()) * F[loc].ns());
239
240// vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > PsiLran(N_sites);
241// for(size_t i=0; i<N_sites; i++)
242// {
243// PsiLran[i] = (F[i].ns() * F[i].c(G[i]));
244// }
245// vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > PsiRran(N_sites);
246// for(size_t i=0; i<N_sites; i++)
247// {
248// PsiRran[i] = (F[i].c(G[i]) * F[i].ns());
249// }
250//
251// vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > first {PsidagLloc,PsidagRloc};
252// vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > last {PsiRran,PsiLran};
253// push_full("Xfull", "Xᵢⱼ", first, last, {-std::sqrt(2.) * std::sqrt(2.), -std::sqrt(2.) * std::sqrt(2.)}, PROP::FERMIONIC);
254// }
255
256 // Local terms: Hubbard-U, t⟂, V⟂, J⟂
257
258 param1d Uph = P.fill_array1d<double>("Uph", "Uphorb", orbitals, loc%Lcell);
259 param2d tperp = P.fill_array2d<double>("tRung", "t", "tPerp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
260 param2d Vperp = P.fill_array2d<double>("Vrung", "V", "Vperp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
261 param2d Jperp = P.fill_array2d<double>("Jrung", "J", "Jperp", orbitals, loc%Lcell, P.get<bool>("CYLINDER"));
262
263 labellist[loc].push_back(Uph.label);
264 labellist[loc].push_back(tperp.label);
265 labellist[loc].push_back(Vperp.label);
266 labellist[loc].push_back(Jperp.label);
267
268 auto Hloc = Mpo<Symmetry_,double>::get_N_site_interaction(F[loc].HubbardHamiltonian(Uph.a, tperp.a, Vperp.a, Jperp.a));
269 pushlist.push_back(std::make_tuple(loc, Hloc, 1.));
270
271 // Nearest-neighbour terms: t, V, J, X
272
273// if (!P.HAS("tFull"))
274// {
275// param2d tpara = P.fill_array2d<double>("t", "tPara", {orbitals, next_orbitals}, loc%Lcell);
276// labellist[loc].push_back(tpara.label);
277//
278// if (loc < N_sites-1 or !static_cast<bool>(boundary))
279// {
280// for (std::size_t alfa=0; alfa<orbitals; ++alfa)
281// for (std::size_t beta=0; beta<next_orbitals; ++beta)
282// {
283// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> cdag_sign_local = (F[loc].cdag(G[loc], alfa) * F[loc].sign());
284// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> c_tight = F[lp1].c(G[lp1], beta);
285
286// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(cdag_sign_local, c_tight), -std::sqrt(2.)*std::sqrt(2.)*tpara(alfa,beta)));
287// }
288// }
289// }
290//
291// if (!P.HAS("Vfull"))
292// {
293// param2d Vpara = P.fill_array2d<double>("V", "Vpara", {orbitals, next_orbitals}, loc%Lcell);
294// labellist[loc].push_back(Vpara.label);
295//
296// if (loc < N_sites-1 or !static_cast<bool>(boundary))
297// {
298// for (std::size_t alfa=0; alfa<orbitals; ++alfa)
299// for (std::size_t beta=0; beta<next_orbitals; ++beta)
300// {
301// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> Tdag_local = F[loc].Tdag(alfa);
302// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> T_tight = F[lp1].T(beta);
303//
304// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(Tdag_local, T_tight), std::sqrt(3.)*Vpara(alfa,beta)));
305// }
306// }
307// }
308//
309// if (!P.HAS("Jfull"))
310// {
311// param2d Jpara = P.fill_array2d<double>("J", "Jpara", {orbitals, next_orbitals}, loc%Lcell);
312// labellist[loc].push_back(Jpara.label);
313//
314// if (loc < N_sites-1 or !static_cast<bool>(boundary))
315// {
316// for (std::size_t alfa=0; alfa<orbitals; ++alfa)
317// for (std::size_t beta=0; beta<next_orbitals; ++beta)
318// {
319// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> Sdag_local = F[loc].Sdag(alfa);
320// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> S_tight = F[lp1].S(beta);
321
322// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(Sdag_local, S_tight), std::sqrt(3.)*Jpara(alfa,beta)));
323// }
324// }
325// }
326//
327// if (!P.HAS("Xfull"))
328// {
329// param2d Xpara = P.fill_array2d<double>("X", "Xpara", {orbitals, next_orbitals}, loc%Lcell);
330// labellist[loc].push_back(Xpara.label);
331//
332// if (loc < N_sites-1 or !static_cast<bool>(boundary))
333// {
334// for (std::size_t alfa=0; alfa<orbitals; ++alfa)
335// for (std::size_t beta=0; beta<next_orbitals; ++beta)
336// {
337// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiRdag_loc = ((F[loc].ns(alfa) * F[loc].cdag(G[loc],alfa)) * F[loc].sign());
338// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiR_tight = (F[lp1].c(G[lp1],beta) * F[lp1].ns(beta));
339//
340// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiLdag_loc = ((F[loc].cdag(G[loc],alfa) * F[loc].sign()) * F[loc].ns(alfa));
341// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiL_tight = (F[lp1].ns(beta) * F[lp1].c(G[lp1],beta));
342
343// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(PsiLdag_loc, PsiR_tight), -std::sqrt(2.)*std::sqrt(2.)*Xpara(alfa,beta)));
344// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(PsiRdag_loc, PsiL_tight), -std::sqrt(2.)*std::sqrt(2.)*Xpara(alfa,beta)));
345// }
346// }
347// }
348//
349// if (!P.HAS("tFull"))
350// {
351// // tPrimePrime
352// param2d tPrimePrime = P.fill_array2d<double>("tPrimePrime", "tPrimePrime_array", {orbitals, next3_orbitals}, loc%Lcell);
353// labellist[loc].push_back(tPrimePrime.label);
354//
355// if (loc < N_sites-2 or !static_cast<bool>(boundary))
356// {
357// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> sign_tight = F[lp1].sign();
358// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> sign_nextn = F[lp2].sign();
359//
360// for (std::size_t alfa=0; alfa<orbitals; ++alfa)
361// for (std::size_t beta=0; beta<next3_orbitals; ++beta)
362// {
363// SiteOperatorQ<Symmetry_, Eigen::MatrixXd> cdag_sign_local = (F[loc].cdag(G[loc],alfa) * F[loc].sign());
364// SiteOperatorQ<Symmetry_, Eigen::MatrixXd> c_nnextn = F[lp3].c(G[lp3],beta);
365//
366// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(cdag_sign_local, sign_tight, sign_nextn, c_nnextn), -std::sqrt(2.)*std::sqrt(2.)*tPrimePrime(alfa,beta)));
367// }
368// }
369// }
370 }
371}
372
373} //end namespace VMPS
374#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_CHARGE_SU2(), Mpo< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::SU2< Sym::ChargeSU2 > >, double > >::type T(size_t locx, size_t locy=0, double factor=1.) const
std::enable_if< Dummy::IS_SPIN_SU2() and!Dummy::IS_CHARGE_SU2(), Mpo< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::SU2< Sym::ChargeSU2 > >, double > >::type P(size_t locx1, size_t locx2, size_t locy1=0, size_t locy2=0) const
Mpo< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::SU2< Sym::ChargeSU2 > >, double > Id() const
vector< FermionBase< Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::SU2< Sym::ChargeSU2 > > > > F
std::size_t N_phys
Definition: MpoTerms.h:400
std::size_t N_sites
Definition: MpoTerms.h:395
DMRG::VERBOSITY::OPTION VERB
Definition: MpoTerms.h:102
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 precalc_TwoSiteData(bool FORCE=false)
bool HAMILTONIAN
Definition: Mpo.h:160
bool HERMITIAN
Definition: Mpo.h:159
Definition: SU2.h:36
Eigen::Index Index
Definition: HubbardU1xSU2.h:26
static void set_operators(const std::vector< FermionBase< Symmetry_ > > &F, const vector< SUB_LATTICE > &G, const ParamHandler &P, PushType< SiteOperator< Symmetry_, double >, double > &pushlist, std::vector< std::vector< std::string > > &labellist, const BC boundary=BC::OPEN)
static constexpr MODEL_FAMILY FAMILY
Definition: HubbardU1xSU2.h:75
HubbardU1xSU2(Mpo< Symmetry > &Mpo_input, const vector< Param > &params)
Definition: HubbardU1xSU2.h:37
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixType
Definition: HubbardU1xSU2.h:28
Symmetry::qType qType
Definition: HubbardU1xSU2.h:27
Sym::S1xS2< Sym::U1< Sym::SpinU1 >, Sym::SU2< Sym::ChargeSU2 > > Symmetry
Definition: HubbardU1xSU2.h:21
static const map< string, any > sweep_defaults
Definition: HubbardU1xSU2.h:79
static const map< string, any > defaults
Definition: HubbardU1xSU2.h:78
static constexpr int spinfac
Definition: HubbardU1xSU2.h:76
static qarray< 2 > singlet(int N=0, int L=0)
Definition: HubbardU1xSU2.h:69
SiteOperatorQ< Symmetry, MatrixType > OperatorType
Definition: HubbardU1xSU2.h:30
#define MAKE_TYPEDEFS(MODEL)
Definition: macros.h:4
const bool COMPRESS
Definition: DmrgTypedefs.h:499
const bool NON_UNITARY
Definition: DmrgTypedefs.h:495
const bool FERMIONIC
Definition: DmrgTypedefs.h:496
const bool HERMITIAN
Definition: DmrgTypedefs.h:492
Definition: qarray.h:26