1#ifndef HUBBARDMODELSU2CHARGE_H_
2#define HUBBARDMODELSU2CHARGE_H_
10#include "Geometry2D.h"
47 typedef Eigen::Matrix<
double,Eigen::Dynamic,Eigen::Dynamic>
MatrixType;
66 size_t Lcell =
P.size();
68 for (
size_t l=0; l<
N_sites; ++l)
N_phys +=
P.get<
size_t>(
"Ly",l%Lcell);
69 this->
calc(P.get<
size_t>(
"maxPower"));
78 template<
typename Symmetry_>
81 const BC boundary=BC::OPEN);
95 {
"t",1.}, {
"tPrime",0.}, {
"tRung",1.}, {
"tPrimePrime",0.},
98 {
"V",0.}, {
"Vrung",0.},
99 {
"Jz",0.}, {
"Jzrung",0.}, {
"Jxy",0.}, {
"Jxyrung",0.},
100 {
"J",0.}, {
"Jperp",0.},
101 {
"Bz",0.}, {
"Bx",0.},
102 {
"X",0.}, {
"Xrung",0.},
103 {
"REMOVE_DOUBLE",
false}, {
"REMOVE_EMPTY",
false}, {
"REMOVE_UP",
false}, {
"REMOVE_DN",
false}, {
"mfactor",1}, {
"k",0},
104 {
"maxPower",2ul}, {
"CYLINDER",
false}, {
"Ly",1ul}
109 {
"max_alpha",100.}, {
"min_alpha",1.}, {
"lim_alpha",11ul}, {
"eps_svd",1e-7},
110 {
"Dincr_abs", 4ul}, {
"Dincr_per", 2ul}, {
"Dincr_rel", 1.1},
111 {
"min_Nsv",0ul}, {
"max_Nrich",-1},
112 {
"max_halfsweeps",24ul}, {
"min_halfsweeps",6ul},
113 {
"Minit",1ul}, {
"Qinit",1ul}, {
"Mlimit",500ul},
114 {
"tol_eigval",1e-7}, {
"tol_state",1e-6},
124 ParamHandler P(params,defaults);
125 size_t Lcell = P.size();
127 for (
size_t l=0; l<N_sites; ++l)
129 N_phys += P.get<
size_t>(
"Ly",l%Lcell);
130 setLocBasis(F[l].get_basis().qloc(),l);
133 param1d U = P.fill_array1d<
double>(
"U",
"Uorb", F[0].orbitals(), 0);
134 if (isfinite(U.a.sum()))
136 this->set_name(
"Hubbard");
140 this->set_name(
"U=∞-Hubbard");
144 std::vector<std::vector<std::string>> labellist;
148 this->construct_from_pushlist(pushlist, labellist, Lcell);
151 this->precalc_TwoSiteData();
154template<
typename Symmetry_>
158 std::size_t Lcell =
P.size();
162 for (std::size_t loc=0; loc<
N_sites; ++loc)
168 std::size_t orbitals =
F[loc].orbitals();
169 std::size_t next_orbitals =
F[lp1].orbitals();
170 std::size_t nextn_orbitals =
F[lp2].orbitals();
171 std::size_t nnextn_orbitals =
F[lp3].orbitals();
187 ss <<
"Ly=" <<
P.get<
size_t>(
"Ly",loc%Lcell);
188 labellist[loc].push_back(ss.str());
190 auto push_full = [&
N_sites, &loc, &
F, &
P, &pushlist, &labellist, &boundary] (
string xxxFull,
string label,
191 const vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > &first,
192 const vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > &last,
193 vector<double> factor,
bool FERMIONIC) ->
void
195 ArrayXXd Full =
P.get<Eigen::ArrayXXd>(xxxFull);
196 vector<vector<std::pair<size_t,double> > > R = Geometry2D::rangeFormat(Full);
198 if (
static_cast<bool>(boundary)) {assert(R.size() ==
N_sites and
"Use an (N_sites)x(N_sites) hopping matrix for open BC!");}
199 else {assert(R.size() >= 2*
N_sites and
"Use at least a (2*N_sites)x(N_sites) hopping matrix for infinite BC!");}
201 for (
size_t j=0; j<first.size(); j++)
202 for (
size_t h=0; h<R[loc].size(); ++h)
204 size_t range = R[loc][h].first;
205 double value = R[loc][h].second;
209 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > ops(range+1);
211 for (
size_t i=1; i<range; ++i)
213 if (FERMIONIC) {ops[i] =
F[(loc+i)%
N_sites].sign();}
216 ops[range] = last[j][(loc+range)%
N_sites];
217 pushlist.push_back(std::make_tuple(loc, ops, factor[j] * value));
222 ss <<
label <<
"(" << Geometry2D::hoppingInfo(Full) <<
")";
223 labellist[loc].push_back(ss.str());
229 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > cup_ranges(
N_sites);
231 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > cdn_ranges(
N_sites);
232 for (
size_t i=0; i<
N_sites; i++)
235 cup_ranges[i] =
F[i].c(
UP,
G[i],0);
237 for (
size_t i=0; i<
N_sites; i++)
240 cdn_ranges[i] =
F[i].c(
DN,
G[i],0);
243 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > frst {cdagup_sign_local, cdagdn_sign_local};
244 vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > last {cup_ranges, cdn_ranges};
245 push_full(
"tFull",
"tᵢⱼ", frst, last, {-std::sqrt(2.),-std::sqrt(2.)},
PROP::FERMIONIC);
249 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > first {
F[loc].Sp(0),
F[loc].Sm(0),
F[loc].Sz(0)};
250 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > Sp_ranges(
N_sites);
251 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > Sm_ranges(
N_sites);
252 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > Sz_ranges(
N_sites);
253 for (
size_t i=0; i<
N_sites; i++)
255 Sp_ranges[i] =
F[i].Sp(0);
256 Sm_ranges[i] =
F[i].Sm(0);
257 Sz_ranges[i] =
F[i].Sz(0);
260 vector<vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > > last {Sm_ranges, Sp_ranges, Sz_ranges};
261 push_full(
"Jfull",
"Jᵢⱼ", first, last, {0.5,0.5,1.},
PROP::BOSONIC);
265 param1d Uph =
P.fill_array1d<
double>(
"Uph",
"Uphorb", orbitals, loc%Lcell);
266 param1d V =
P.fill_array1d<
double>(
"V",
"Vorb", orbitals, loc%Lcell);
269 param2d tPerp =
P.fill_array2d<
double>(
"tRung",
"t",
"tPerp", orbitals, loc%Lcell,
P.get<
bool>(
"CYLINDER"));
270 param2d Jxyperp =
P.fill_array2d<
double>(
"JxyRung",
"Jxy",
"JxyPerp", orbitals, loc%Lcell,
P.get<
bool>(
"CYLINDER"));
271 param2d Jzperp =
P.fill_array2d<
double>(
"JzRung",
"Jz",
"JzPerp", orbitals, loc%Lcell,
P.get<
bool>(
"CYLINDER"));
272 param2d Jperp =
P.fill_array2d<
double>(
"JRung",
"J",
"JPerp", orbitals, loc%Lcell,
P.get<
bool>(
"CYLINDER"));
273 param1d Bz =
P.fill_array1d<
double>(
"Bz",
"Bzorb", orbitals, loc%Lcell);
274 param1d Bx =
P.fill_array1d<
double>(
"Bx",
"Bxorb", orbitals, loc%Lcell);
276 labellist[loc].push_back(Uph.label);
279 labellist[loc].push_back(tPerp.label);
280 labellist[loc].push_back(Jxyperp.label);
281 labellist[loc].push_back(Jzperp.label);
282 labellist[loc].push_back(Jperp.label);
283 labellist[loc].push_back(Bz.label);
284 labellist[loc].push_back(Bx.label);
286 ArrayXXd Vperp =
F[loc].ZeroHopping();
288 auto sum_array = [] (
const ArrayXXd& a1,
const ArrayXXd& a2)
290 ArrayXXd res(a1.rows(), a1.cols());
291 for (
int i=0; i<a1.rows(); ++i)
292 for (
int j=0; j<a1.rows(); ++j)
294 res(i,j) = a1(i,j) + a2(i,j);
302 F[loc].
template HubbardHamiltonian<double>(Uph.a, tPerp.a, Vperp, sum_array(Jperp.a,Jzperp.a), sum_array(Jperp.a,Jxyperp.a), Bz.a, Bx.a)
304 pushlist.push_back(std::make_tuple(loc, Hloc, 1.));
std::enable_if< Dummy::IS_SPIN_SU2() and!Dummy::IS_CHARGE_SU2(), Mpo< Sym::SU2< Sym::ChargeSU2 >, double > >::type P(size_t locx1, size_t locx2, size_t locy1=0, size_t locy2=0) const
Mpo< Sym::SU2< Sym::ChargeSU2 >, double > Id() const
vector< FermionBase< Sym::SU2< Sym::ChargeSU2 > > > F
DMRG::VERBOSITY::OPTION VERB
void calc(const std::size_t power)
static std::vector< T > get_N_site_interaction(T const &Op0, Operator const &... Ops)
void precalc_TwoSiteData(bool FORCE=false)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixType
HubbardSU2charge(Mpo< Symmetry > &Mpo_input, const vector< Param > ¶ms)
static constexpr MODEL_FAMILY FAMILY
Sym::SU2< Sym::ChargeSU2 > Symmetry
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 const map< string, any > sweep_defaults
static constexpr int spinfac
static qarray< 1 > singlet(int N=0)
static const map< string, any > defaults
#define MAKE_TYPEDEFS(MODEL)
void finalize(bool PRINT_STATS=false)