1#ifndef HUBBARDMODELU1XU1_H_COMPLEX
2#define HUBBARDMODELU1XU1_H_COMPLEX
10#include "Geometry2D.h"
15 public HubbardObservables<Sym::S1xS2<Sym::U1<Sym::SpinU1>,Sym::U1<Sym::ChargeU1> >,complex<double>>,
22 typedef Eigen::Matrix<complex<
double>,Eigen::Dynamic,Eigen::Dynamic>
MatrixType;
40 size_t Lcell =
P.size();
42 for (
size_t l=0; l<
N_sites; ++l)
N_phys +=
P.get<
size_t>(
"Ly",l%Lcell);
59 template<
typename Symmetry_>
61 PushType<
SiteOperator<Symmetry_,complex<double>>,complex<double>>& pushlist, std::vector<std::vector<std::string>>& labellist,
62 const BC boundary=BC::OPEN);
76 {
"t",1.+0.i}, {
"tPrime",0.i}, {
"tRung",1.+0.i}, {
"tPrimePrime",0.i},
79 {
"V",0.}, {
"Vext",0.}, {
"Vrung",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}
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},
106 size_t Lcell =
P.size();
108 for (
size_t l=0; l<
N_sites; ++l)
110 N_phys +=
P.get<
size_t>(
"Ly",l%Lcell);
117 std::vector<std::vector<std::string>> labellist;
126template<
typename Symmetry_>
130 std::size_t Lcell =
P.size();
134 for (std::size_t loc=0; loc<
N_sites; ++loc)
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();
146 ss <<
"Ly=" <<
P.get<
size_t>(
"Ly",loc%Lcell);
147 labellist[loc].push_back(ss.str());
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,
154 bool FERMIONIC) ->
void
156 ArrayXXcd Full =
P.get<Eigen::ArrayXXcd>(xxxFull);
157 vector<vector<std::pair<size_t,complex<double>> > > R = Geometry2D::rangeFormat(Full);
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!");}
162 for (
size_t j=0; j<first.size(); j++)
163 for (
size_t h=0; h<R[loc].size(); ++h)
165 size_t range = R[loc][h].first;
166 complex<double> value = R[loc][h].second;
170 vector<SiteOperatorQ<Symmetry_,MatrixType> > ops(range+1);
172 for (
size_t i=1; i<range; ++i)
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>>();}
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));
185 ss <<
label <<
"(" << Geometry2D::hoppingInfo(Full) <<
")";
186 labellist[loc].push_back(ss.str());
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>>();}
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);
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++)
212 Tz_ranges[i] =
F[i].Tz(0).template cast<complex<double>>();
215 vector<vector<SiteOperatorQ<Symmetry_,MatrixType> > > last {Tz_ranges};
216 push_full(
"Vzfull",
"Vzᵢⱼ", first, last, {1.}, {
false},
PROP::BOSONIC);
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++)
225 tz_ranges[i] =
F[i].tz(0).template cast<complex<double>>();
228 vector<vector<SiteOperatorQ<Symmetry_,MatrixType> > > last {tz_ranges};
229 push_full(
"vzfull",
"vzᵢⱼ", first, last, {1.}, {
false},
PROP::BOSONIC);
232 if (
P.HAS(
"VextFull"))
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);
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++)
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>>();
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);
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"));
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();
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>>(),
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>>())
297 pushlist.push_back(std::make_tuple(loc, Hloc, 1.+0.i));
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
void finalize(const bool COMPRESS=true, const std::size_t power=1, const double tolerance=::mynumeric_limits< double >::epsilon())
void setLocBasis(const std::vector< std::vector< qType > > &q)
DMRG::VERBOSITY::OPTION VERB
MpoTerms< Symmetry, OtherScalar > cast()
void set_name(const std::string &label_in)
static std::vector< T > get_N_site_interaction(T const &Op0, Operator const &... Ops)
void precalc_TwoSiteData(bool FORCE=false)
void construct_from_pushlist(const PushType< OperatorType, CouplScalar > &pushlist, const std::vector< std::vector< std::string > > &labellist, size_t Lcell)
static constexpr int spinfac
static constexpr MODEL_FAMILY FAMILY
PeierlsHubbardU1xU1(Mpo< Symmetry, complex< double > > &Mpo_input, const vector< Param > ¶ms)
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)