VMPS++
Loading...
Searching...
No Matches
HubbardSU2charge.h
Go to the documentation of this file.
1#ifndef HUBBARDMODELSU2CHARGE_H_
2#define HUBBARDMODELSU2CHARGE_H_
3
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
39class HubbardSU2charge : public Mpo<Sym::SU2<Sym::ChargeSU2> ,double>,
40 public HubbardObservables<Sym::SU2<Sym::ChargeSU2> >,
41 public ParamReturner
42{
43public:
44
47 typedef Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic> MatrixType;
49
50private:
51
52 typedef Eigen::Index Index;
54
55public:
56
59
60 HubbardSU2charge(Mpo<Symmetry> &Mpo_input, const vector<Param> &params)
61 :Mpo<Symmetry>(Mpo_input),
64 {
65 ParamHandler P(params,HubbardSU2charge::defaults);
66 size_t Lcell = P.size();
67 N_phys = 0;
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"));
70 this->precalc_TwoSiteData();
71 this->HERMITIAN = true;
72 this->HAMILTONIAN = true;
73 };
74
75 HubbardSU2charge (const size_t &L, const vector<Param> &params, const BC &boundary=BC::OPEN, const DMRG::VERBOSITY::OPTION &VERB=DMRG::VERBOSITY::OPTION::ON_EXIT);
77
78 template<typename Symmetry_>
79 static void set_operators (const std::vector<FermionBase<Symmetry_> > &F, const vector<SUB_LATTICE> &G, const ParamHandler &P,
80 PushType<SiteOperator<Symmetry_,double>,double>& pushlist, std::vector<std::vector<std::string>>& labellist,
81 const BC boundary=BC::OPEN);
82
83 static qarray<1> singlet (int N=0) {return qarray<1>{1};};
84 static constexpr MODEL_FAMILY FAMILY = HUBBARD;
85 static constexpr int spinfac = 1;
86
87 static const map<string,any> defaults;
88 static const map<string,any> sweep_defaults;
89};
90
91// V is standard next-nearest neighbour density interaction
92// Vz and Vxy are anisotropic isospin-isospin next-nearest neighbour interaction
93const map<string,any> HubbardSU2charge::defaults =
94{
95 {"t",1.}, {"tPrime",0.}, {"tRung",1.}, {"tPrimePrime",0.},
96 {"mu",0.}, {"t0",0.},
97 {"U",0.}, {"Uph",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}
105};
106
107const map<string,any> HubbardSU2charge::sweep_defaults =
108{
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},
115 {"savePeriod",0ul}, {"CALC_S_ON_EXIT", true}, {"CONVTEST", DMRG::CONVTEST::VAR_2SITE}
116};
117
119HubbardSU2charge (const size_t &L, const vector<Param> &params, const BC &boundary, const DMRG::VERBOSITY::OPTION &VERB)
120:Mpo<Symmetry> (L, qarray<Symmetry::Nq>({1}), "", PROP::HERMITIAN, PROP::NON_UNITARY, boundary, VERB),
123{
124 ParamHandler P(params,defaults);
125 size_t Lcell = P.size();
126
127 for (size_t l=0; l<N_sites; ++l)
128 {
129 N_phys += P.get<size_t>("Ly",l%Lcell);
130 setLocBasis(F[l].get_basis().qloc(),l);
131 }
132
133 param1d U = P.fill_array1d<double>("U", "Uorb", F[0].orbitals(), 0);
134 if (isfinite(U.a.sum()))
135 {
136 this->set_name("Hubbard");
137 }
138 else
139 {
140 this->set_name("U=∞-Hubbard");
141 }
142
144 std::vector<std::vector<std::string>> labellist;
145 HubbardSU2charge::set_operators(F, G, P, pushlist, labellist, boundary); // F, G are set in HubbardObservables
146 //add_operators(F, P, pushlist, labellist, boundary);
147
148 this->construct_from_pushlist(pushlist, labellist, Lcell);
149 this->finalize(PROP::COMPRESS, P.get<size_t>("maxPower"));
150
151 this->precalc_TwoSiteData();
152}
153
154template<typename Symmetry_>
156set_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)
157{
158 std::size_t Lcell = P.size();
159 std::size_t N_sites = F.size();
160 if(labellist.size() != N_sites) {labellist.resize(N_sites);}
161
162 for (std::size_t loc=0; loc<N_sites; ++loc)
163 {
164 size_t lp1 = (loc+1)%N_sites;
165 size_t lp2 = (loc+2)%N_sites;
166 size_t lp3 = (loc+3)%N_sites;
167
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();
172
173// vector<SUB_LATTICE> G(N_sites);
174// if (P.HAS("G")) {G = P.get<vector<SUB_LATTICE> >("G");}
175// else // set default (-1)^l
176// {
177// G[0] = static_cast<SUB_LATTICE>(1);
178// for (int l=1; l<N_sites; l+=1) G[l] = static_cast<SUB_LATTICE>(-1*G[l-1]);
179// }
180
181// auto Gloc = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,loc)));
182// auto Glp1 = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,lp1)));
183// auto Glp2 = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,lp2)));
184// auto Glp3 = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,lp3)));
185
186 stringstream ss;
187 ss << "Ly=" << P.get<size_t>("Ly",loc%Lcell);
188 labellist[loc].push_back(ss.str());
189
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
194 {
195 ArrayXXd Full = P.get<Eigen::ArrayXXd>(xxxFull);
196 vector<vector<std::pair<size_t,double> > > R = Geometry2D::rangeFormat(Full);
197
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!");}
200
201 for (size_t j=0; j<first.size(); j++)
202 for (size_t h=0; h<R[loc].size(); ++h)
203 {
204 size_t range = R[loc][h].first;
205 double value = R[loc][h].second;
206
207 if (range != 0)
208 {
209 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > ops(range+1);
210 ops[0] = first[j];
211 for (size_t i=1; i<range; ++i)
212 {
213 if (FERMIONIC) {ops[i] = F[(loc+i)%N_sites].sign();}
214 else {ops[i] = F[(loc+i)%N_sites].Id();}
215 }
216 ops[range] = last[j][(loc+range)%N_sites];
217 pushlist.push_back(std::make_tuple(loc, ops, factor[j] * value));
218 }
219 }
220
221 stringstream ss;
222 ss << label << "(" << Geometry2D::hoppingInfo(Full) << ")";
223 labellist[loc].push_back(ss.str());
224 };
225
226 if (P.HAS("tFull"))
227 {
228 SiteOperatorQ<Symmetry_,Eigen::MatrixXd> cdagup_sign_local = F[loc].cdag(UP,G[loc],0) * F[loc].sign();
229 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > cup_ranges(N_sites);
230 SiteOperatorQ<Symmetry_,Eigen::MatrixXd> cdagdn_sign_local = F[loc].cdag(DN,G[loc],0) * F[loc].sign();
231 vector<SiteOperatorQ<Symmetry_,Eigen::MatrixXd> > cdn_ranges(N_sites);
232 for (size_t i=0; i<N_sites; i++)
233 {
234 //auto Gi = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,i)));
235 cup_ranges[i] = F[i].c(UP,G[i],0);
236 }
237 for (size_t i=0; i<N_sites; i++)
238 {
239 //auto Gi = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,i)));
240 cdn_ranges[i] = F[i].c(DN,G[i],0);
241 }
242
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);
246 }
247 if (P.HAS("Jfull"))
248 {
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++)
254 {
255 Sp_ranges[i] = F[i].Sp(0);
256 Sm_ranges[i] = F[i].Sm(0);
257 Sz_ranges[i] = F[i].Sz(0);
258 }
259
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);
262 }
263 // Local terms: U, t0, μ, t⟂, V⟂, J⟂
264
265 param1d Uph = P.fill_array1d<double>("Uph", "Uphorb", orbitals, loc%Lcell);
266 param1d V = P.fill_array1d<double>("V", "Vorb", orbitals, loc%Lcell);
267 //param1d t0 = P.fill_array1d<double>("t0", "t0orb", orbitals, loc%Lcell);
268 //param1d mu = P.fill_array1d<double>("mu", "muorb", 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);
275
276 labellist[loc].push_back(Uph.label);
277 //labellist[loc].push_back(t0.label);
278 //labellist[loc].push_back(mu.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);
285
286 ArrayXXd Vperp = F[loc].ZeroHopping();
287
288 auto sum_array = [] (const ArrayXXd& a1, const ArrayXXd& a2)
289 {
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)
293 {
294 res(i,j) = a1(i,j) + a2(i,j);
295 }
296 return res;
297 };
298
300 (
301 //HubbardHamiltonian(U.a,tPerp.a,Vperp,Jzsubperp,Jxysubperp,Bzsub.a,Bxsub.a));
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)
303 );
304 pushlist.push_back(std::make_tuple(loc, Hloc, 1.));
305
306 // Nearest-neighbour terms: t, V, J
307
308// if (!P.HAS("tFull") and !P.HAS("Jfull"))
309// {
310// param2d tpara = P.fill_array2d<double>("t", "tPara", {orbitals, next_orbitals}, loc%Lcell);
311// param2d Vpara = P.fill_array2d<double>("V", "Vpara", {orbitals, next_orbitals}, loc%Lcell);
312// param2d Vzpara = P.fill_array2d<double>("Vz", "Vzpara", {orbitals, next_orbitals}, loc%Lcell);
313// param2d Vxypara = P.fill_array2d<double>("Vxy", "Vxypara", {orbitals, next_orbitals}, loc%Lcell);
314// param2d Jpara = P.fill_array2d<double>("J", "Jpara", {orbitals, next_orbitals}, loc%Lcell);
315// param2d Xpara = P.fill_array2d<double>("X", "Xpara", {orbitals, next_orbitals}, loc%Lcell);
316//
317// labellist[loc].push_back(tpara.label);
318// labellist[loc].push_back(Vpara.label);
319// labellist[loc].push_back(Vzpara.label);
320// labellist[loc].push_back(Vxypara.label);
321// labellist[loc].push_back(Jpara.label);
322// labellist[loc].push_back(Xpara.label);
323//
324// if (loc < N_sites-1 or !static_cast<bool>(boundary))
325// {
326// for (std::size_t alfa=0; alfa<orbitals; ++alfa)
327// for (std::size_t beta=0; beta<next_orbitals; ++beta)
328// {
329// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> c_sign_local = (F[loc].c(alfa) * F[loc].sign());
330// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> cdag_sign_local = (F[loc].cdag(alfa) * F[loc].sign());
331//
332// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> c_tight = F[lp1].c (beta);
333// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> cdag_tight = F[lp1].cdag(beta);
334//
335// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> n_local = F[loc].n(alfa);
336// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> n_tight = F[lp1].n(beta);
337//
338// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> tz_local = F[loc].Tz(alfa);
339// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> tz_tight = F[lp1].Tz(beta);
340//
341// auto Gloc = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,loc)));
342// auto Glp1 = static_cast<SUB_LATTICE>(static_cast<int>(pow(-1,lp1)));
343// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> tp_local = F[loc].Tp(alfa,Gloc);
344// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> tm_tight = F[lp1].Tm(beta,Glp1);
345//
346// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> tm_local = F[loc].Tm(alfa,Gloc);
347// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> tp_tight = F[lp1].Tp(beta,Glp1);
348//
349// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> Sdag_local = F[loc].Sdag(alfa);
350// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> S_tight = F[lp1].S (beta);
351//
352// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiLloc = ((F[loc].ns(alfa) * F[loc].c(alfa)) * F[loc].sign());
353// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiRloc = ((F[loc].c(alfa) * F[loc].sign()) * F[loc].ns(alfa));
354// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiLlp1 = (F[lp1].ns(beta) * F[lp1].c(beta));
355// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsiRlp1 = (F[lp1].c(beta) * F[lp1].ns(beta));
356// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsidagLloc = ((F[loc].cdag(alfa) * F[loc].sign()) * F[loc].ns(alfa));
357// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsidagRloc = ((F[loc].ns(alfa) * F[loc].cdag(alfa)) * F[loc].sign());
358// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsidagLlp1 = (F[lp1].cdag(beta) * F[lp1].ns(beta));
359// SiteOperatorQ<Symmetry_,Eigen::MatrixXd> PsidagRlp1 = (F[lp1].ns(beta) * F[lp1].cdag(beta));
360//
361// //hopping
362// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(cdag_sign_local, c_tight), -std::sqrt(2.)*tpara(alfa,beta)));
363// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(c_sign_local, cdag_tight), -std::sqrt(2.)*tpara(alfa,beta)));
364//
365// //density-density interaction
366// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(n_local, n_tight), Vpara(alfa,beta)));
367//
368// //isospin-isopsin interaction
369// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(tp_local, tm_tight), 0.5*Vxypara(alfa,beta)));
370// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(tm_local, tp_tight), 0.5*Vxypara(alfa,beta)));
371// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(tz_local, tz_tight), Vzpara (alfa,beta)));
372//
373// //spin-spin interaction
374// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(Sdag_local, S_tight), std::sqrt(3.)*Jpara(alfa,beta)));
375//
376// //correlated hopping
377// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(PsidagLloc, PsiRlp1), -std::sqrt(2.)*Xpara(alfa,beta)));
378// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(PsidagRloc, PsiLlp1), -std::sqrt(2.)*Xpara(alfa,beta)));
379// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(PsiLloc, PsidagRlp1), -std::sqrt(2.)*Xpara(alfa,beta)));
380// pushlist.push_back(std::make_tuple(loc, Mpo<Symmetry_,double>::get_N_site_interaction(PsiRloc, PsidagLlp1), -std::sqrt(2.)*Xpara(alfa,beta)));
381// }
382// }
383// }
384 }
385}
386
387} // end namespace VMPS::models
388
389#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::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
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
void calc(const std::size_t power)
Definition: MpoTerms.h:1135
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::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixType
HubbardSU2charge(Mpo< Symmetry > &Mpo_input, const vector< Param > &params)
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)
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
void finalize(bool PRINT_STATS=false)
Definition: functions.h:127
Definition: qarray.h:26