ABACUS/include/JSC_Heis.h

778 lines
38 KiB
C++

/**********************************************************
This software is part of J.-S. Caux's ABACUS++ library.
Copyright (c)
-----------------------------------------------------------
File: Heis.h
Purpose: Declares Heisenberg chain classes and functions.
***********************************************************/
#ifndef _HEIS_
#define _HEIS_
#include "JSC.h"
namespace JSC {
// First, some global constants...
const long long int ID_UPPER_LIMIT = 10000000LL; // max size of vectors we can define without seg fault
const int INTERVALS_SIZE = 100000; // size of Scan_Intervals arrays
const int NBASESMAX = 1000; // max number of bases kept
const DP ITER_REQ_PREC = 100.0 * MACHINE_EPS_SQ;
//const DP ITER_REQ_PREC = MACHINE_EPS_SQ;
// Cutoffs on particle numbers
//const int NPARTICLES_MAX = 24;
//const int NHOLES_MAX = NPARTICLES_MAX/2;
//const int MAX_RAPS_ABOVE_ZERO = 10; // max nr of rapidities above lowest type
//const int NPARTICLES_MAX = 2;
//const int NHOLES_MAX = 1;
//const int MAX_TYPES_IN_BASE = 4; // maximal number of particle types we allow in bases
const int MAXSTRINGS = 20; // maximal number of particle types we allow in bases
//const int MAXSTRINGS = 2; // maximal number of particle types we allow in bases
//const DP HEIS_deltaprec = 1.0e-6;//1.0e-4; // maximal string deviation allowed // DEPRECATED in ++T_9
const int NEXC_MAX_HEIS = 16; // maximal number of excitations (string binding/unbinding, particle-hole) considered
//***********************************************************************
class Heis_Chain {
public:
DP J;
DP Delta;
DP anis; // acos(Delta) if Delta < 1.0, 0 if Delta == 1.0, acosh(Delta) if Delta > 1.0
DP hz;
int Nsites;
int Nstrings; // how many possible strings. The following two arrays have Nstrings nonzero elements.
int* Str_L; // vector (length M) containing the allowed string lengths. Elements that are 0 have no meaning.
int* par; // vector (length M) containing the parities of the strings. Elements that are 0 have no meaning.
// Parities are all +1 except for gapless XXZ subcases
DP* si_n_anis_over_2; // for optimization: sin for XXZ, sinh for XXZ_gpd
DP* co_n_anis_over_2; // for optimization
DP* ta_n_anis_over_2; // for optimization
DP prec; // precision required for computations, always put to ITER_REQ_PREC
public:
Heis_Chain ();
Heis_Chain (DP JJ, DP DD, DP hh, int NN); // contructor: simply initializes
Heis_Chain (const Heis_Chain& RefChain); // copy constructor;
Heis_Chain& operator= (const Heis_Chain& RefChain);
bool operator== (const Heis_Chain& RefChain);
bool operator!= (const Heis_Chain& RefChain);
~Heis_Chain(); // destructor
public:
//void Scan_for_Possible_Bases (int Mdown, Vect<long long int>& possible_base_id, int& nfound, int nexc_max_used,
// int base_level_to_scan, Vect<int>& Nrapidities);
//Vect<long long int> Possible_Bases (int Mdown); // returns a vector of possible bases
/* Deactivated in ++G_8
void Scan_for_Possible_Bases (int Mdown, Vect<string>& possible_base_label, int& nfound, int nexc_max_used,
int base_level_to_scan, Vect<int>& Nrapidities);
Vect<string> Possible_Bases (int Mdown); // returns a vector of possible bases
*/
};
//****************************************************************************
// Objects in class Heis_Base are a checked vector containing the number of rapidities of allowable types for a given state
class Heis_Base {
public:
int Mdown; // total number of down spins
Vect<int> Nrap; // Nrap[i] contains the number of rapidities of type i, i = 0, Nstrings - 1.
int Nraptot; // total number of strings in this state
Vect<DP> Ix2_infty; // Ix2_infty[i] contains the max of BAE function for the (half-)integer I[i], i = 0, Nstrings - 1.
Vect<int> Ix2_min;
Vect<int> Ix2_max; // Ix2_max[i] contains the integer part of 2*I_infty, with correct parity for base.
//long long int id; // identification number
double dimH; // dimension of sub Hilbert space associated to this base; use double to avoid max int problems.
string baselabel; // base label
public:
Heis_Base ();
Heis_Base (const Heis_Base& RefBase); // copy constructor
Heis_Base (const Heis_Chain& RefChain, int M); // constructs configuration with all Mdown in one-string of +1 parity
Heis_Base (const Heis_Chain& RefChain, const Vect<int>& Nrapidities); // sets to Nrapidities vector, and checks consistency
//Heis_Base (const Heis_Chain& RefChain, long long int id_ref);
Heis_Base (const Heis_Chain& RefChain, string baselabel_ref);
inline int& operator[] (const int i);
inline const int& operator[] (const int i) const;
Heis_Base& operator= (const Heis_Base& RefBase);
bool operator== (const Heis_Base& RefBase);
bool operator!= (const Heis_Base& RefBase);
void Compute_Ix2_limits(const Heis_Chain& RefChain); // computes the Ix2_infty and Ix2_max
//void Scan_for_Possible_Types (Vect<long long int>& possible_type_id, int& nfound, int base_level, Vect<int>& Nexcitations);
//Vect<long long int> Possible_Types (); // returns a vector of possible types
};
inline int& Heis_Base::operator[] (const int i)
{
return Nrap[i];
}
inline const int& Heis_Base::operator[] (const int i) const
{
return Nrap[i];
}
//****************************************************************************
/*
// Objects in class Ix2_Config carry all the I's of a given state
class Ix2_Config {
//private:
public:
int Nstrings;
Vect<int> Nrap;
int Nraptot;
//int** Ix2;
Vect<Vect<int> > Ix2;
public:
Ix2_Config ();
Ix2_Config (const Heis_Chain& RefChain, int M); // constructor, puts I's to ground state
Ix2_Config (const Heis_Chain& RefChain, const Heis_Base& base); // constructor, putting I's to lowest-energy config
// consistent with Heis_Base configuration for chain RefChain
Ix2_Config& operator= (const Ix2_Config& RefConfig);
//inline int* operator[] (const int i);
inline Vect<int> operator[] (const int i);
//inline const int* operator[] (const int i) const;
inline const Vect<int> operator[] (const int i) const;
//~Ix2_Config(); // not needed, inherited from Vect
string Return_Label (const Ix2_Config& OriginIx2);
};
//inline int* Ix2_Config::operator[] (const int i)
inline Vect<int> Ix2_Config::operator[] (const int i)
{
return Ix2[i];
}
//inline const int* Ix2_Config::operator[] (const int i) const
inline const Vect<int> Ix2_Config::operator[] (const int i) const
{
return Ix2[i];
}
std::ostream& operator<< (std::ostream& s, const Ix2_Config& RefConfig);
*/
//****************************************************************************
// Objects in class Lambda carry all rapidities of a state
class Lambda {
private:
int Nstrings;
Vect<int> Nrap;
int Nraptot;
DP** lambda;
//Vect<Vect<DP> > lambda;
public:
Lambda ();
Lambda (const Heis_Chain& RefChain, int M); // constructor, puts all lambda's to zero
Lambda (const Heis_Chain& RefChain, const Heis_Base& base); // constructor, putting I's to lowest-energy config
// consistent with Heis_Base configuration for chain RefChain
Lambda& operator= (const Lambda& RefConfig);
inline DP* operator[] (const int i);
//inline Vect<DP> operator[] (const int i);
inline const DP* operator[] (const int i) const;
//inline const Vect<DP> operator[] (const int i) const;
~Lambda();
};
inline DP* Lambda::operator[] (const int i)
//inline Vect<DP> Lambda::operator[] (const int i)
{
return lambda[i];
}
inline const DP* Lambda::operator[] (const int i) const
//inline const Vect<DP> Lambda::operator[] (const int i) const
{
return lambda[i];
}
//****************************************************************************
// Objects in class Ix2_Offsets carry Young tableau representations of the Ix2 configurations
/*
class Ix2_Offsets {
public:
Heis_Base base;
Vect<Young_Tableau> Tableau; // vector of pointers to tableaux at each level
long long int type_id;
long long int id; // id number of offset
long long int maxid; // max id number allowable
public:
Ix2_Offsets ();
Ix2_Offsets (const Ix2_Offsets& RefOffset); // copy constructor
Ix2_Offsets (const Heis_Base& RefBase, long long int req_type_id);
Ix2_Offsets (const Heis_Base& RefBase, Vect<int> nparticles); // sets all tableaux to empty ones, with nparticles[] at each level
Ix2_Offsets& operator= (const Ix2_Offsets& RefOffset);
bool operator<= (const Ix2_Offsets& RefOffsets);
bool operator>= (const Ix2_Offsets& RefOffsets);
public:
void Set_to_id (long long int idnr);
void Compute_id ();
void Compute_type_id ();
public:
bool Add_Boxes_From_Lowest (int Nboxes, bool odd_sectors); // adds Nboxes in minimal energy config, all boxes in either even or odd sectors
};
inline long long int Ix2_Offsets_type_id (Vect<int>& nparticles)
{
long long int type_id_here = 0ULL;
for (int i = 0; i < nparticles.size(); ++i)
type_id_here += nparticles[i] * pow_ulli(10ULL, i);
return(type_id_here);
}
long long int Find_idmin (Ix2_Offsets& scan_offsets, int particle_type, int tableau_level, int Row_L_min);
long long int Find_idmax (Ix2_Offsets& scan_offsets, int particle_type, int tableau_level);
*/
//****************************************************************************
// Objects in class Ix2_Offsets_List carry a vector of used Ix2_Offsets
/*
class Ix2_Offsets_List {
public:
int ndef;
Vect<Ix2_Offsets> Offsets;
public:
Ix2_Offsets_List ();
Ix2_Offsets& Return_Offsets (Heis_Base& RefBase, Vect<int> nparticles); // returns the Ix2_Offsets corresponding to nparticles[]/base
Ix2_Offsets& Return_Offsets (Heis_Base& RefBase, long long int req_type_id);
};
*/
//****************************************************************************
// Objects in class Heis_Bethe_State carry all information about an eigenstate
// Derived classes include XXZ_Bethe_State, XXX_Bethe_State, XXZ_gpd_Bethe_State
// These contain subclass-specific functions and data.
class Heis_Bethe_State {
public:
Heis_Chain chain;
Heis_Base base;
//Ix2_Offsets offsets;
//Ix2_Config Ix2;
Vect<Vect<int> > Ix2;
Lambda lambda;
Lambda deviation; // string deviations
Lambda BE; // Bethe equation for relevant rapidity, in the form BE = theta - (1/N)\sum ... - \pi I/N = 0
DP diffsq; // sum of squares of rapidity differences in last iteration
int conv; // convergence status
DP dev; // sum of absolute values of string deviations
int iter; // number of iterations necessary for convergence
int iter_Newton; // number of iterations necessary for convergence (Newton method)
DP E; // total energy
int iK; // K = 2.0*PI * iK/Nsites
DP K; // total momentum
DP lnnorm; // ln of norm of reduced Gaudin matrix
//long long int base_id;
//long long int type_id;
//long long int id;
//long long int maxid;
//int nparticles;
string label;
public:
Heis_Bethe_State ();
Heis_Bethe_State (const Heis_Bethe_State& RefState); // copy constructor
//Heis_Bethe_State (const Heis_Bethe_State& RefState, long long int type_id_ref); // new state with requested type_id
Heis_Bethe_State (const Heis_Chain& RefChain, int M); // constructor to ground-state configuration
Heis_Bethe_State (const Heis_Chain& RefChain, const Heis_Base& base); // constructor to lowest-energy config with base
//Heis_Bethe_State (const Heis_Chain& RefChain, long long int base_id_ref, long long int type_id_ref);
virtual ~Heis_Bethe_State () {};
public:
int Charge () { return(base.Mdown); };
//void Set_Ix2_Offsets (const Ix2_Offsets& RefOffset); // sets the Ix2 to given offsets
//void Set_to_id (long long int id_ref);
//void Set_to_id (long long int id_ref, Heis_Bethe_State& RefState);
//int Nparticles (); // counts the number of particles in state once Ix2 offsets set (so type_id is correctly set)
//void Set_to_Label (string label_ref, const Ix2_Config& OriginIx2);
void Set_to_Label (string label_ref, const Vect<Vect<int> >& OriginIx2);
void Set_Label_from_Ix2 (const Vect<Vect<int> >& OriginIx2);
bool Check_Symmetry (); // checks whether the I's are symmetrically distributed
void Compute_diffsq (); // \sum BE[j][alpha]^2
void Find_Rapidities (bool reset_rapidities); // Finds the rapidities
void Find_Rapidities_Twisted (bool reset_rapidities, DP twist); // Finds the rapidities with twist added to RHS of logBE
//void BAE_smackdown (DP max_allowed);
//void Solve_BAE_smackdown (DP max_allowed, int maxruns);
void Solve_BAE_bisect (int j, int alpha, DP req_prec, int itermax);
void Iterate_BAE (DP iter_factor); // Finds new set of lambda[j][alpha] from previous one by simple iteration
void Solve_BAE_straight_iter (DP straight_prec, int max_iter_interp, DP iter_factor);
void Solve_BAE_extrap (DP extrap_prec, int max_iter_extrap, DP iter_factor);
void Iterate_BAE_Newton (); // Finds new set of lambda[j][alpha] from previous one by a Newton step
void Solve_BAE_Newton (DP Newton_prec, int max_iter_Newton);
void Solve_BAE_with_silk_gloves (DP silk_prec, int max_iter_silk, DP iter_factor);
void Compute_lnnorm ();
void Compute_Momentum ();
void Compute_All (bool reset_rapidities); // solves BAE, computes E, K and lnnorm
inline bool Set_to_Inner_Skeleton (int iKneeded, const Vect<Vect<int> >& OriginStateIx2)
{
Ix2[0][0] = Ix2[0][1] - 2;
Ix2[0][base.Nrap[0] - 1] = Ix2[0][base.Nrap[0] - 2] + 2;
(*this).Compute_Momentum();
if (base.Nrap[0] == 0) return(false);
if (iKneeded >= iK) Ix2[0][base.Nrap[0]-1] += 2*(iKneeded - iK);
else Ix2[0][0] += 2*(iKneeded - iK);
if (Ix2[0][0] < base.Ix2_min[0] || Ix2[0][base.Nrap[0]-1] > base.Ix2_max[0]) return(false);
(*this).Set_Label_from_Ix2 (OriginStateIx2);
return(true);
}
void Set_to_Outer_Skeleton (const Vect<Vect<int> >& OriginStateIx2) {
Ix2[0][0] = base.Ix2_min[0] - 4;
Ix2[0][base.Nrap[0]-1] = base.Ix2_max[0] + 4;
(*this).Set_Label_from_Ix2 (OriginStateIx2);
};
void Set_to_Closest_Matching_Ix2_fixed_Base (const Heis_Bethe_State& StateToMatch); // defined in Heis.cc
// Virtual functions, all defined in the derived classes
public:
virtual void Set_Free_lambdas() { JSCerror("Heis_Bethe_State::..."); } // sets the rapidities to solutions of BAEs without scattering terms
virtual bool Check_Admissibility(char option) { JSCerror("Heis_Bethe_State::..."); return(false); }
// verifies that we don't have a symmetrical Ix2 config with a Ix2 == 0 for a string of even length >= 2.
virtual void Compute_BE (int j, int alpha) { JSCerror("Heis_Bethe_State::..."); }
virtual void Compute_BE () { JSCerror("Heis_Bethe_State::..."); }
virtual DP Iterate_BAE(int i, int alpha) { JSCerror("Heis_Bethe_State::..."); return(0.0);}
virtual bool Check_Rapidities() { JSCerror("Heis_Bethe_State::..."); return(false); }
virtual DP String_delta () { JSCerror("Heis_Bethe_State::..."); return(0.0); }
virtual void Compute_Energy () { JSCerror("Heis_Bethe_State::..."); }
virtual void Build_Reduced_Gaudin_Matrix (SQMat<complex<DP> >& Gaudin_Red) { JSCerror("Heis_Bethe_State::..."); }
};
inline bool Is_Inner_Skeleton (Heis_Bethe_State& State) {
return (State.base.Nrap[0] >= 2 && (State.Ix2[0][0] == State.Ix2[0][1] - 2 || State.Ix2[0][State.base.Nrap[0]-1] == State.Ix2[0][State.base.Nrap[0]-2] + 2));
};
inline bool Is_Outer_Skeleton (Heis_Bethe_State& State) {
return (State.Ix2[0][0] == State.base.Ix2_min[0] - 4 && State.Ix2[0][State.base.Nrap[0]-1] == State.base.Ix2_max[0] + 4);
};
inline bool Force_Descent (char whichDSF, Heis_Bethe_State& ScanState, Heis_Bethe_State& RefState, int desc_type_required, int iKmod, DP Chem_Pot)
{
bool force_descent = false;
// Force descent if energy of ScanState is lower than that of RefState
//if (ScanState.E - RefState.E - (ScanState.base.Mdown - RefState.base.Mdown) < 0.0) return(true);
/*
// We force descent if
// 1) - there exists a higher string whose quantum number is still on 0
// AND - there is at most a single particle-hole in the 0 base level
// AND - either the particle or the hole hasn't yet moved.
if (ScanState.base_id/100000LL > 0) { // there is a higher string
int type0 = ScanState.type_id % 10000;
if (type0 == 0
|| type0 == 101 && ScanState.offsets.Tableau[0].id * ScanState.offsets.Tableau[2].id == 0LL
|| type0 == 110 && ScanState.offsets.Tableau[1].id * ScanState.offsets.Tableau[2].id == 0LL
|| type0 == 1001 && ScanState.offsets.Tableau[0].id * ScanState.offsets.Tableau[3].id == 0LL
|| type0 == 1010 && ScanState.offsets.Tableau[1].id * ScanState.offsets.Tableau[3].id == 0LL) // single p-h pair in base level 0
for (int j = 1; j < ScanState.chain.Nstrings; ++j) {
if (ScanState.base[j] == 1 && ScanState.Ix2[j][0] == 0) {
force_descent = true;
}
}
}
*/
// Force descent if quantum nr distribution is symmetric:
//if (ScanState.Check_Symmetry()) force_descent = true;
//if (desc_type_required > 8 && ScanState.Check_Symmetry()) force_descent = true;
// Force descent for longitudinal if we're at zero or pi momentum:
//ScanState.Compute_Momentum();
//if (whichDSF == 'z' && (ScanState.iK - RefState.iK) % iKmod == 0) force_descent = true;
//if (desc_type_required > 8 && whichDSF == 'z' && (2*(ScanState.iK - RefState.iK) % iKmod == 0)) force_descent = true; // type_req > 8 means that we don't conserve momentum
// Force descent for all DSFs if we're at K = 0 or PI and not conserving momentum upon descent:
if (desc_type_required > 8 && (2*(ScanState.iK - RefState.iK) % iKmod == 0)) force_descent = true; // type_req > 8 means that we don't conserve momentum
//if (force_descent) cout << "Forcing descent on state with label " << ScanState.label << endl;
if (ScanState.chain.Delta > 1.0) {
/*
// Count the nr of holes in one-strings:
int nholes = 0;
for (int i = 0; i < ScanState.base.Nrap[0] - 1; ++i) if (ScanState.Ix2[0][i+1] - ScanState.Ix2[0][i] != 2) nholes++;
if (nholes <= 2) {
if (ScanState.base.Nrap[0] == ScanState.base.Mdown - 2 && ScanState.base.Nrap[1] == 1) force_descent = true;
if (ScanState.base.Nrap[0] == ScanState.base.Mdown - 3 && ScanState.base.Nrap[2] == 1) force_descent = true;
if (ScanState.base.Nrap[0] == ScanState.base.Mdown - 4 && ScanState.base.Nrap[1] == 2) force_descent = true;
}
*/
if (ScanState.label.compare(0, 10, "14x1y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "14x1y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "14x1y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "14x1y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "14x1y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "14x1y1_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_0x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_1x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "12x1y2_2x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "13x2y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "13x2y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "13x2y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "13x2y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "13x2y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "13x2y1_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "30x1y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "30x1y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "30x1y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "30x1y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "30x1y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "30x1y1_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_0x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_1x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "28x1y2_2x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "29x2y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "29x2y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "29x2y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "29x2y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "29x2y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "29x2y1_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "46x1y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "46x1y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "46x1y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "46x1y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "46x1y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "46x1y1_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_0x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_1x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "44x1y2_2x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "45x2y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "45x2y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "45x2y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "45x2y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "45x2y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "45x2y1_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "62x1y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "62x1y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "62x1y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "62x1y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "62x1y1_2x0") == 0) force_descent = true;
//if (ScanState.label.compare(0, 10, "62x1y1_2x1") == 0 && desc_type_required < 2) force_descent = true;
if (ScanState.label.compare(0, 10, "62x1y1_2x1") == 0
&& (desc_type_required == 14 || desc_type_required == 13 || desc_type_required == 11 || desc_type_required == 10)) force_descent = true;
/*
if (ScanState.label.compare(0, 10, "60x1y2_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_0x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_1x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_2x1") == 0 && desc_type_required < 2) force_descent = true;
if (ScanState.label.compare(0, 10, "60x1y2_2x2") == 0 && desc_type_required < 2) force_descent = true;
if (ScanState.label.compare(0, 10, "61x2y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "61x2y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "61x2y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "61x2y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "61x2y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "61x2y1_2x1") == 0 && desc_type_required < 2) force_descent = true;
*/
if (ScanState.label.compare(0, 11, "126x1y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "126x1y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "126x1y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "126x1y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "126x1y1_2x0") == 0) force_descent = true;
//if (ScanState.label.compare(0, 11, "126x1y1_2x1") == 0 && desc_type_required < 2) force_descent = true;
if (ScanState.label.compare(0, 11, "126x1y1_2x1") == 0
&& (desc_type_required == 14 || desc_type_required == 13 || desc_type_required == 11 || desc_type_required == 10)) force_descent = true;
/*
if (ScanState.label.compare(0, 10, "124x1y2_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_0x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_1x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_2x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "124x1y2_2x2") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "125x2y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "125x2y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "125x2y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "125x2y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "125x2y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 10, "125x2y1_2x1") == 0) force_descent = true;
*/
if (ScanState.label.compare(0, 11, "254x1y1_0x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "254x1y1_0x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "254x1y1_1x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "254x1y1_1x1") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "254x1y1_2x0") == 0) force_descent = true;
if (ScanState.label.compare(0, 11, "254x1y1_2x1") == 0 && desc_type_required < 2) force_descent = true;
// Do not force descent a state with rapidities outside of fundamental interval:
/*
for (int j = 0; j < ScanState.chain.Nstrings; ++j) {
// Logic: allow rapidities -PI/2 <= lambda <= PI/2 (including boundaries)
if (ScanState.base.Nrap[j] > 0 &&
(ScanState.lambda[j][0] < -PI/2 || ScanState.lambda[j][ScanState.base.Nrap[j] - 1] > PI/2))
force_descent = false;
// rapidities should also be ordered as the quantum numbers:
for (int alpha = 1; alpha < ScanState.base.Nrap[j]; ++alpha)
if (ScanState.lambda[j][alpha - 1] >= ScanState.lambda[j][alpha])
force_descent = false;
}
*/
//if (force_descent) cout << "Forcing descent on state with label " << ScanState.label << endl;
} // if Delta > 1
//if (ScanState.base.Nrap[0] == ScanState.base.Mdown - 2 && ScanState.base.Nrap[1] == 1 && ScanState.Ix2[1][0] == 0) force_descent = true;
//if (ScanState.base.Nrap[0] == ScanState.base.Mdown - 3 && ScanState.base.Nrap[2] == 1 && ScanState.Ix2[2][0] == 0) force_descent = true;
//if (ScanState.base.Nrap[0] == ScanState.base.Mdown - 4 && ScanState.base.Nrap[1] == 2 && ScanState.Ix2[1][0] == -1 && ScanState.Ix2[1][1] == 1) force_descent = true;
return(force_descent);
}
std::ostream& operator<< (std::ostream& s, const Heis_Bethe_State& state);
//****************************************************************************
// Objects in class XXZ_Bethe_State carry all extra information pertaining to XXZ gapless
class XXZ_Bethe_State : public Heis_Bethe_State {
public:
Lambda sinhlambda;
Lambda coshlambda;
Lambda tanhlambda;
public:
XXZ_Bethe_State ();
XXZ_Bethe_State (const XXZ_Bethe_State& RefState); // copy constructor
XXZ_Bethe_State (const Heis_Chain& RefChain, int M); // constructor to ground-state configuration
XXZ_Bethe_State (const Heis_Chain& RefChain, const Heis_Base& base); // constructor to lowest-energy config with base
//XXZ_Bethe_State (const Heis_Chain& RefChain, long long int base_id_ref, long long int type_id_ref); // constructor to lowest-energy config with bas
public:
XXZ_Bethe_State& operator= (const XXZ_Bethe_State& RefState);
public:
void Set_Free_lambdas(); // sets the rapidities to solutions of BAEs without scattering terms
void Compute_sinhlambda();
void Compute_coshlambda();
void Compute_tanhlambda();
bool Check_Admissibility(char option); // verifies that we don't have a symmetrical Ix2 config with a Ix2 == 0 for a string of even length >= 2.
void Compute_BE (int j, int alpha);
void Compute_BE ();
DP Iterate_BAE(int i, int j);
bool Check_Rapidities(); // checks that all rapidities are not nan
DP String_delta ();
void Compute_Energy ();
//void Compute_Momentum ();
void Build_Reduced_Gaudin_Matrix (SQMat<complex<DP> >& Gaudin_Red);
// XXZ specific functions:
public:
};
XXZ_Bethe_State Add_Particle_at_Center (const XXZ_Bethe_State& RefState);
XXZ_Bethe_State Remove_Particle_at_Center (const XXZ_Bethe_State& RefState);
//****************************************************************************
// Objects in class XXX_Bethe_State carry all extra information pertaining to XXX antiferromagnet
class XXX_Bethe_State : public Heis_Bethe_State {
public:
XXX_Bethe_State ();
XXX_Bethe_State (const XXX_Bethe_State& RefState); // copy constructor
XXX_Bethe_State (const Heis_Chain& RefChain, int M); // constructor to ground-state configuration
XXX_Bethe_State (const Heis_Chain& RefChain, const Heis_Base& base); // constructor to lowest-energy config with base
//XXX_Bethe_State (const Heis_Chain& RefChain, long long int base_id_ref, long long int type_id_ref); // constructor to lowest-energy config with base
public:
XXX_Bethe_State& operator= (const XXX_Bethe_State& RefState);
public:
void Set_Free_lambdas(); // sets the rapidities to solutions of BAEs without scattering terms
bool Check_Admissibility(char option); // verifies that we don't have a symmetrical Ix2 config with a Ix2 == 0 for a string of even length >= 2.
void Compute_BE (int j, int alpha);
void Compute_BE ();
DP Iterate_BAE(int i, int j);
bool Check_Rapidities(); // checks that all rapidities are not nan
DP String_delta ();
void Compute_Energy ();
//void Compute_Momentum ();
void Build_Reduced_Gaudin_Matrix (SQMat<complex<DP> >& Gaudin_Red);
// XXX specific functions
public:
bool Check_Finite_rap ();
};
XXX_Bethe_State Add_Particle_at_Center (const XXX_Bethe_State& RefState);
XXX_Bethe_State Remove_Particle_at_Center (const XXX_Bethe_State& RefState);
//****************************************************************************
// Objects in class XXZ_gpd_Bethe_State carry all extra information pertaining to XXZ gapped antiferromagnets
class XXZ_gpd_Bethe_State : public Heis_Bethe_State {
public:
Lambda sinlambda;
Lambda coslambda;
Lambda tanlambda;
public:
XXZ_gpd_Bethe_State ();
XXZ_gpd_Bethe_State (const XXZ_gpd_Bethe_State& RefState); // copy constructor
XXZ_gpd_Bethe_State (const Heis_Chain& RefChain, int M); // constructor to ground-state configuration
XXZ_gpd_Bethe_State (const Heis_Chain& RefChain, const Heis_Base& base); // constructor to lowest-energy config with base
//XXZ_gpd_Bethe_State (const Heis_Chain& RefChain, long long int base_id_ref, long long int type_id_ref); // constructor to lowest-energy config with base
public:
XXZ_gpd_Bethe_State& operator= (const XXZ_gpd_Bethe_State& RefState);
public:
void Set_Free_lambdas(); // sets the rapidities to solutions of BAEs without scattering terms
void Compute_sinlambda();
void Compute_coslambda();
void Compute_tanlambda();
int Weight(); // weight function for contributions cutoff
bool Check_Admissibility(char option); // verifies that we don't have a symmetrical Ix2 config with a Ix2 == 0 for a string of even length >= 2.
void Compute_BE (int j, int alpha);
void Compute_BE ();
DP Iterate_BAE(int i, int j);
void Iterate_BAE_Newton();
bool Check_Rapidities(); // checks that all rapidities are not nan and are in interval ]-PI/2, PI/2]
DP String_delta ();
void Compute_Energy ();
//void Compute_Momentum ();
void Build_Reduced_Gaudin_Matrix (SQMat<complex<DP> >& Gaudin_Red);
// XXZ_gpd specific functions
public:
};
XXZ_gpd_Bethe_State Add_Particle_at_Center (const XXZ_gpd_Bethe_State& RefState);
XXZ_gpd_Bethe_State Remove_Particle_at_Center (const XXZ_gpd_Bethe_State& RefState);
//***********************************************
// Function declarations
// in M_vs_H.cc
DP Ezero (DP Delta, int N, int M);
DP H_vs_M (DP Delta, int N, int M);
DP HZmin (DP Delta, int N, int M, Vect_DP& Ezero_ref);
int M_vs_H (DP Delta, int N, DP HZ);
DP X_avg (char xyorz, DP Delta, int N, int M);
DP Chemical_Potential (const Heis_Bethe_State& RefState);
DP Particle_Hole_Excitation_Cost (char whichDSF, Heis_Bethe_State& AveragingState);
//DP Sumrule_Factor (char whichDSF, Heis_Bethe_State& RefState, DP Chem_Pot, bool fixed_iK, int iKneeded);
DP Sumrule_Factor (char whichDSF, Heis_Bethe_State& RefState, DP Chem_Pot, int iKmin, int iKmax);
void Evaluate_F_Sumrule (string prefix, char whichDSF, const Heis_Bethe_State& RefState, DP Chem_Pot, int iKmin, int iKmax);
complex<DP> ln_Sz_ME (XXZ_Bethe_State& A, XXZ_Bethe_State& B);
complex<DP> ln_Smin_ME (XXZ_Bethe_State& A, XXZ_Bethe_State& B);
complex<DP> ln_Sz_ME (XXX_Bethe_State& A, XXX_Bethe_State& B);
complex<DP> ln_Smin_ME (XXX_Bethe_State& A, XXX_Bethe_State& B);
// From Antoine Klauser:
complex<DP> ln_Szz_ME (XXX_Bethe_State& A, XXX_Bethe_State& B);
complex<DP> ln_Szm_p_Smz_ME (XXX_Bethe_State& A, XXX_Bethe_State& B);
complex<DP> ln_Smm_ME (XXX_Bethe_State& A, XXX_Bethe_State& B);
complex<DP> ln_Sz_ME (XXZ_gpd_Bethe_State& A, XXZ_gpd_Bethe_State& B);
complex<DP> ln_Smin_ME (XXZ_gpd_Bethe_State& A, XXZ_gpd_Bethe_State& B);
// The following functions have become member functions.
//DP String_delta (XXZ_Bethe_State& state);
//DP String_delta (XXX_Bethe_State& state);
//DP String_delta (XXZ_gpd_Bethe_State& state);
//DP Compute_Form_Factor_Entry (char whichDSF, Heis_Bethe_State& LeftState, Heis_Bethe_State& RefState, DP Chem_Pot, fstream& DAT_outfile);
//DP Compute_Matrix_Element_Contrib (char whichDSF, int iKmin, int iKmax, XXZ_Bethe_State& LeftState,
// XXZ_Bethe_State& RefState, DP Chem_Pot, fstream& DAT_outfile);
DP Compute_Matrix_Element_Contrib (char whichDSF, int iKmin, int iKmax, XXZ_Bethe_State& LeftState,
XXZ_Bethe_State& RefState, DP Chem_Pot, stringstream& DAT_outfile);
//DP Compute_Matrix_Element_Contrib (char whichDSF, int iKmin, int iKmax, XXX_Bethe_State& LeftState,
// XXX_Bethe_State& RefState, DP Chem_Pot, fstream& DAT_outfile);
DP Compute_Matrix_Element_Contrib (char whichDSF, int iKmin, int iKmax, XXX_Bethe_State& LeftState,
XXX_Bethe_State& RefState, DP Chem_Pot, stringstream& DAT_outfile);
//DP Compute_Matrix_Element_Contrib (char whichDSF, int iKmin, int iKmax, XXZ_gpd_Bethe_State& LeftState,
// XXZ_gpd_Bethe_State& RefState, DP Chem_Pot, fstream& DAT_outfile);
DP Compute_Matrix_Element_Contrib (char whichDSF, int iKmin, int iKmax, XXZ_gpd_Bethe_State& LeftState,
XXZ_gpd_Bethe_State& RefState, DP Chem_Pot, stringstream& DAT_outfile);
// For geometric quench:
complex<DP> ln_Overlap (XXX_Bethe_State& A, XXX_Bethe_State& B);
void Scan_Heis_Geometric_Quench (DP Delta, int N_1, int M, long long int base_id_1, long long int type_id_1, long long int id_1,
int N_2, int iKmin, int iKmax, int Max_Secs, bool refine);
} // namespace JSC
#endif