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