ABACUS/src/HEIS/XXZ_gpd_Bethe_State.cc

815 строки
28 KiB
C++

/**********************************************************
This software is part of J.-S. Caux's ABACUS++ library.
Copyright (c)
-----------------------------------------------------------
File: src/HEIS/XXZ_gpd_Bethe_State.cc
Purpose: Defines all functions for XXZ_gpd_Bethe_State
******************************************************************/
#include "JSC.h"
using namespace std;
namespace JSC {
// Function prototypes
inline DP fbar_XXZ_gpd (DP lambda, int par, DP tanhnetaover2);
DP Theta_XXZ_gpd (DP lambda, int nj, int nk, DP* tanhnetaover2);
DP hbar_XXZ_gpd (DP lambda, int n, DP* si_n_anis_over_2);
DP ddlambda_Theta_XXZ_gpd (DP lambda, int nj, int nk, DP* si_n_anis_over_2);
//***************************************************************************************************
// Function definitions: class XXZ_gpd_Bethe_State
XXZ_gpd_Bethe_State::XXZ_gpd_Bethe_State ()
: Heis_Bethe_State(), sinlambda(Lambda(chain, 1)), coslambda(Lambda(chain, 1)), tanlambda(Lambda(chain, 1))
{};
XXZ_gpd_Bethe_State::XXZ_gpd_Bethe_State (const XXZ_gpd_Bethe_State& RefState) // copy constructor
: Heis_Bethe_State(RefState),
sinlambda(Lambda(RefState.chain, RefState.base)), coslambda(Lambda(RefState.chain, RefState.base)),
tanlambda(Lambda(RefState.chain, RefState.base))
{
// copy arrays into new ones
for (int j = 0; j < RefState.chain.Nstrings; ++j) {
for (int alpha = 0; alpha < RefState.base[j]; ++j) {
sinlambda[j][alpha] = RefState.sinlambda[j][alpha];
coslambda[j][alpha] = RefState.coslambda[j][alpha];
tanlambda[j][alpha] = RefState.tanlambda[j][alpha];
}
}
}
XXZ_gpd_Bethe_State::XXZ_gpd_Bethe_State (const Heis_Chain& RefChain, int M)
: Heis_Bethe_State(RefChain, M),
sinlambda(Lambda(RefChain, M)), coslambda(Lambda(RefChain, M)), tanlambda(Lambda(RefChain, M))
{
if (RefChain.Delta <= 1.0) JSCerror("Delta too low in XXZ_gpd_Bethe_State constructor");
}
XXZ_gpd_Bethe_State::XXZ_gpd_Bethe_State (const Heis_Chain& RefChain, const Heis_Base& RefBase)
: Heis_Bethe_State(RefChain, RefBase),
sinlambda(Lambda(RefChain, RefBase)), coslambda(Lambda(RefChain, RefBase)),
tanlambda(Lambda(RefChain, RefBase))
{
if (RefChain.Delta <= 1.0) JSCerror("Delta too low in XXZ_gpd_Bethe_State constructor");
}
/*
XXZ_gpd_Bethe_State::XXZ_gpd_Bethe_State (const Heis_Chain& RefChain, long long int base_id_ref, long long int type_id_ref)
: Heis_Bethe_State(RefChain, base_id_ref, type_id_ref),
sinlambda(Lambda(chain, base)), coslambda(Lambda(chain, base)), tanlambda(Lambda(chain, base))
{
if (RefChain.Delta <= 1.0) JSCerror("Delta too low in XXZ_gpd_Bethe_State constructor");
}
*/
XXZ_gpd_Bethe_State& XXZ_gpd_Bethe_State::operator= (const XXZ_gpd_Bethe_State& RefState)
{
if (this != &RefState) {
chain = RefState.chain;
base = RefState.base;
//offsets = RefState.offsets;
Ix2 = RefState.Ix2;
lambda = RefState.lambda;
BE = RefState.BE;
diffsq = RefState.diffsq;
conv = RefState.conv;
iter = RefState.iter;
iter_Newton = RefState.iter_Newton;
E = RefState.E;
iK = RefState.iK;
K = RefState.K;
lnnorm = RefState.lnnorm;
//base_id = RefState.base_id;
//type_id = RefState.type_id;
//id = RefState.id;
//maxid = RefState.maxid;
//nparticles = RefState.nparticles;
label = RefState.label;
sinlambda = RefState.sinlambda;
coslambda = RefState.coslambda;
tanlambda = RefState.tanlambda;
}
return(*this);
}
// Member functions
void XXZ_gpd_Bethe_State::Set_Free_lambdas()
{
// Sets all the rapidities to the solutions of the BAEs without scattering terms
for (int i = 0; i < chain.Nstrings; ++i) {
for (int alpha = 0; alpha < base[i]; ++alpha) {
lambda[i][alpha] = atan((tanh(chain.Str_L[i] * 0.5 * chain.anis) * tan(PI * 0.5 * Ix2[i][alpha]/chain.Nsites)));
}
}
return;
}
void XXZ_gpd_Bethe_State::Compute_sinlambda()
{
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) sinlambda[j][alpha] = sin(lambda[j][alpha]);
}
return;
}
void XXZ_gpd_Bethe_State::Compute_coslambda()
{
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) coslambda[j][alpha] = cos(lambda[j][alpha]);
}
return;
}
void XXZ_gpd_Bethe_State::Compute_tanlambda()
{
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
tanlambda[j][alpha] = tan(lambda[j][alpha]);
//if (lambda[j][alpha] > 0.5*PI) cout << "Rapidity higher than 0.5*PI: j = " << j << "\talpha = " << alpha << "\trap = " << lambda[j][alpha] << "\ttan = " << tanlambda[j][alpha] << endl;
}
}
return;
}
bool XXZ_gpd_Bethe_State::Check_Admissibility(char option)
{
// This function checks the admissibility of the Ix2's of a state:
// returns false if there are higher strings with Ix2 = 0, a totally symmetric distribution of I's at each level,
// and strings of equal length modulo 2 and parity with Ix2 = 0, meaning at least two equal roots in BAE.
bool answer = true;
Vect<bool> Zero_at_level(false, chain.Nstrings); // whether there exists an Ix2 == 0 at a given level
/*
Vect<bool> min_Ix2_max_busy(false, chain.Nstrings);
Vect<bool> plus_Ix2_max_busy(false, chain.Nstrings);
for (int j = 0; j < chain.Nstrings; ++j)
for (int alpha = 0; alpha < base[j]; ++alpha) {
if (Ix2[j][alpha] == -base.Ix2_max[j]) min_Ix2_max_busy[j] = true;
if (Ix2[j][alpha] == base.Ix2_max[j]) plus_Ix2_max_busy[j] = true;
}
*/
/*
// State is not admissible if this is false: -N/2 + 1 \leq \sum I^j_{\alpha} \leq N
int sum_all_Ix2 = 0;
for (int j = 0; j < chain.Nstrings; ++j)
for (int alpha = 0; alpha < base[j]; ++alpha) {
sum_all_Ix2 += Ix2[j][alpha];
}
if (sum_all_Ix2 > chain.Nsites || sum_all_Ix2 <= -chain.Nsites) {
cout << "\tSum Ix2 out of fundamental interval: sum_all_Ix2 = " << sum_all_Ix2 << "\tfor N = " << chain.Nsites << endl;
return(false);
}
*/
//Deactivated 2014 06 11, put in Heis_Form_Factor_Entry.cc
/*
// State is not admissible if I_max > 1/2 (N - \sum_k t_{jk} M_k) or I_min < -1/2 (N - sum_k...) + 1 at any level:
int sum1 = 0;
for (int j = 0; j < chain.Nstrings; ++j) {
sum1 = 0;
for (int k = 0; k < chain.Nstrings; ++k) {
sum1 += base[k] * (2 * JSC::min(chain.Str_L[j], chain.Str_L[k]) - ((j == k) ? 1 : 0));
}
// Define limits...
//if (!((Nrap[j] + Ix2_max[j]) % 2)) Ix2_max[j] -= 1;
// This almost does it: only missing are the states with one on -PI/2 and one on PI/2
if (base[j] >= 1 && (Ix2[j][0] <= -(chain.Nsites - sum1) ||
(Ix2[j][base[j] - 1] - Ix2[j][0]) > 2*(chain.Nsites - sum1))) {
//cout << "\tAn Ix2 is out of interval at level " << j << endl;
//cout << Ix2[j][base[j] - 1] << "\t" << Ix2[j][0] << "\t" << chain.Nsites << "\t" << sum1 << endl;
return(false);
}
}
*/
/*
// State is not admissible if both a min_ and plus_Ix2_max are busy simultaneously:
bool is_a_min_Ix2_max_busy = false;
for (int j = 0; j < chain.Nstrings; ++j) if (min_Ix2_max_busy[j]) is_a_min_Ix2_max_busy = true;
bool is_a_plus_Ix2_max_busy = false;
for (int j = 0; j < chain.Nstrings; ++j) if (plus_Ix2_max_busy[j]) is_a_plus_Ix2_max_busy = true;
*/
/*
// State is not admissible if all min_Ix2_max are busy simultaneously:
bool any_min_Ix2_max_free = false;
for (int j = 0; j < chain.Nstrings; ++j)
if (base[j] > 0 && !min_Ix2_max_busy[j]) any_min_Ix2_max_free = true;
if (!any_min_Ix2_max_free) return(false);
*/
/*
// State is not admissible if -Ix2_max, -Ix2_max + 2, ..., -Ix2_max + 2*(Str_L - 1) are busy:
for (int j = 0; j < chain.Nstrings; ++j)
if (base[j] > 0 && Ix2[j][0] <= -base.Ix2_max[j] + 2*(chain.Str_L[j] - 1))
return(false);
// Almost correct with above !
// State is not admissible if Ix2_max - 2, ..., Ix2_max - 2*(Str_L - 2) are busy (NB: one slot more than on left):
for (int j = 0; j < chain.Nstrings; ++j)
if (base[j] > 0 && Ix2[j][base[j] - 1] >= base.Ix2_max[j] - 2*(chain.Str_L[j] - 2))
return(false);
*/
// Check that at all at least doubly occupied levels, the difference between max and min quantum numbers
// is strictly smaller than 2*Ix2_max - 2, so that lambda_max - lambda_min < PI at each level:
//for (int j = 0; j < chain.Nstrings; ++j)
//if (base[j] >= 2 && Ix2[j][base[j] - 1] - Ix2[j][0] >= 2* base.Ix2_max[j] - 2) answer = false;
//if (is_a_min_Ix2_max_busy && is_a_plus_Ix2_max_busy) return(false);
bool higher_string_on_zero = false;
for (int j = 0; j < chain.Nstrings; ++j) {
// The following line puts answer to true if there is at least one higher string with zero Ix2
for (int alpha = 0; alpha < base[j]; ++alpha) if ((Ix2[j][alpha] == 0) && (chain.Str_L[j] > 2) /*&& !(chain.Str_L[j] % 2)*/)
higher_string_on_zero = true;
for (int alpha = 0; alpha < base[j]; ++alpha) if (Ix2[j][alpha] == 0) Zero_at_level[j] = true;
// NOTE: if base[j] == 0, Zero_at_level[j] remains false.
}
bool symmetric_state = (*this).Check_Symmetry();
bool string_coincidence = false;
// Checks that we have strings of equal length modulo 2 and same parity with Ix2 == 0, so equal rapidities, and inadmissibility
for (int j1 = 0; j1 < chain.Nstrings; ++j1) {
for (int j2 = j1 + 1; j2 < chain.Nstrings; ++j2)
if (Zero_at_level[j1] && Zero_at_level[j2] && (!((chain.Str_L[j1] + chain.Str_L[j2])%2)))
string_coincidence = true;
}
//bool onep_onem_on_zero = false;
//if (option == 'm' || option == 'p') { // for Smin, we also exclude symmetric states with 1+ and 1- strings on zero
//if (Zero_at_level[0] && Zero_at_level[1]) onep_onem_on_zero = true;
//}
//cout << min_Ix2_max_busy << "\t" << symmetric_state << "\t" << higher_string_on_zero << "\t" << string_coincidence << "\t" << onep_onem_on_zero << endl;
//answer = !((symmetric_state && (higher_string_on_zero || string_coincidence || onep_onem_on_zero)));
answer = !(symmetric_state && (higher_string_on_zero || string_coincidence));
// Explicitly exclude state with min_Ix2_max_busy && iK == 0
//Compute_Momentum();
//answer = !((min_Ix2_max_busy && iK == 0) || (symmetric_state && (higher_string_on_zero || string_coincidence || onep_onem_on_zero)));
// Now check that no Ix2 is equal to +N (since we take -N into account, and I + N == I by periodicity of exp)
for (int j = 0; j < chain.Nstrings; ++j)
for (int alpha = 0; alpha < base[j]; ++alpha) if ((Ix2[j][alpha] < -chain.Nsites) || (Ix2[j][alpha] >= chain.Nsites)) answer = false;
if (!answer) {
E = 0.0;
K = 0.0;
conv = 0;
iter = 0;
iter_Newton = 0;
lnnorm = -100.0;
}
return(answer); // answer == true: nothing wrong with this Ix2_config
}
void XXZ_gpd_Bethe_State::Compute_BE (int j, int alpha)
{
// Fills in the BE members with the value of the Bethe equations.
tanlambda[j][alpha] = tan(lambda[j][alpha]);
DP sumtheta = 0.0;
sumtheta = 0.0;
for (int k = 0; k < chain.Nstrings; ++k) {
for (int beta = 0; beta < base[k]; ++beta)
if ((chain.Str_L[j] == 1) && (chain.Str_L[k] == 1)) {
sumtheta += atan ((tanlambda[j][alpha] - tanlambda[k][beta])/((1.0 + tanlambda[j][alpha] * tanlambda[k][beta])
* chain.ta_n_anis_over_2[2]))
+ PI * floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
}
else sumtheta += 0.5 * Theta_XXZ_gpd((tanlambda[j][alpha] - tanlambda[k][beta])/(1.0 + tanlambda[j][alpha] * tanlambda[k][beta]),
chain.Str_L[j], chain.Str_L[k], chain.ta_n_anis_over_2)
+ PI * (2.0 * JSC::min(chain.Str_L[j], chain.Str_L[k]) - ((j == k) ? 1.0 : 0))
* floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
}
sumtheta *= 2.0;
BE[j][alpha] = 2.0 * (atan(tanlambda[j][alpha]/chain.ta_n_anis_over_2[chain.Str_L[j]])
+ PI * floor(0.5 + lambda[j][alpha]/PI))
- (sumtheta + PI*Ix2[j][alpha])/chain.Nsites;
}
void XXZ_gpd_Bethe_State::Compute_BE ()
{
// Fills in the BE members with the value of the Bethe equations.
(*this).Compute_tanlambda();
DP sumtheta = 0.0;
for (int j = 0; j < chain.Nstrings; ++j)
for (int alpha = 0; alpha < base[j]; ++alpha) {
sumtheta = 0.0;
for (int k = 0; k < chain.Nstrings; ++k) {
for (int beta = 0; beta < base[k]; ++beta)
if ((chain.Str_L[j] == 1) && (chain.Str_L[k] == 1)) {
sumtheta += atan ((tanlambda[j][alpha] - tanlambda[k][beta])/((1.0 + tanlambda[j][alpha] * tanlambda[k][beta])
* chain.ta_n_anis_over_2[2]))
+ PI * floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
}
else sumtheta += 0.5 * Theta_XXZ_gpd((tanlambda[j][alpha] - tanlambda[k][beta])/(1.0 + tanlambda[j][alpha] * tanlambda[k][beta]),
chain.Str_L[j], chain.Str_L[k], chain.ta_n_anis_over_2)
+ PI * (2.0 * JSC::min(chain.Str_L[j], chain.Str_L[k]) - ((j == k) ? 1.0 : 0))
* floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
}
sumtheta *= 2.0;
BE[j][alpha] = 2.0 * (atan(tanlambda[j][alpha]/chain.ta_n_anis_over_2[chain.Str_L[j]])
+ PI * floor(0.5 + lambda[j][alpha]/PI))
- (sumtheta + PI*Ix2[j][alpha])/chain.Nsites;
}
}
DP XXZ_gpd_Bethe_State::Iterate_BAE (int j, int alpha)
{
// Returns a new iteration value for lambda[j][alpha] given tanlambda[][] and BE[][]
// Assumes that tanlambda[][] and BE[][] have been computed.
DP arg0 = 0.5 * (2.0 * (atan(tanlambda[j][alpha]/chain.ta_n_anis_over_2[chain.Str_L[j]])
+ PI * floor(0.5 + lambda[j][alpha]/PI)) - BE[j][alpha]);
DP arg = chain.ta_n_anis_over_2[chain.Str_L[j]] * tan(
arg0
//0.5 *
//(PI * Ix2[j][alpha] + sumtheta)/chain.Nsites
//(2.0 * (atan(tanlambda[j][alpha]/chain.ta_n_anis_over_2[chain.Str_L[j]])
// + PI * floor(0.5 + lambda[j][alpha]/PI)) - BE[j][alpha])
);
return(atan(arg)
//+ PI * floor(0.5 + arg0)
//0.5 * (Ix2[j][alpha] + sumtheta/PI)/(chain.Nsites)
+ PI * floor(0.5 + arg0/PI)
);
}
/*
void XXZ_gpd_Bethe_State::Iterate_BAE ()
{
// Recalculates the rapidities by iterating Bethe equations
Lambda New_lambda(chain, base);
DP sumtheta = 0.0;
DP arg = 0.0;
// First, compute the tan of rapidities:
(*this).Compute_tanlambda();
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
sumtheta = 0.0;
for (int k = 0; k < chain.Nstrings; ++k) {
for (int beta = 0; beta < base[k]; ++beta)
if ((chain.Str_L[j] == 1) && (chain.Str_L[k] == 1))
sumtheta += atan((tanlambda[j][alpha] - tanlambda[k][beta])/((1.0 + tanlambda[j][alpha] * tanlambda[k][beta])
* chain.ta_n_anis_over_2[2]))
+ PI * floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
else sumtheta += 0.5 * Theta_XXZ_gpd((tanlambda[j][alpha] - tanlambda[k][beta])/(1.0 + tanlambda[j][alpha] * tanlambda[k][beta]),
chain.Str_L[j], chain.Str_L[k], chain.ta_n_anis_over_2)
+ PI * (2.0 * JSC::min(chain.Str_L[j], chain.Str_L[k]) - ((j == k) ? 1.0 : 0))
* floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
}
sumtheta *= 2.0;
arg = chain.ta_n_anis_over_2[chain.Str_L[j]] * tan((PI * 0.5 * Ix2[j][alpha] + 0.5 * sumtheta)/chain.Nsites);
New_lambda[j][alpha] = atan(arg) + PI * floor(0.5 + (0.5 * Ix2[j][alpha] + 0.5 * sumtheta/PI)/(chain.Nsites));
}
}
DP New_diffsq = 0.0;
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
//New_diffsq += pow(tan(New_lambda[j][alpha]) - tanlambda[j][alpha], 2.0);
New_diffsq += pow(New_lambda[j][alpha] - lambda[j][alpha], 2.0);
lambda[j][alpha] = 1.0 * New_lambda[j][alpha] + 0.0 * lambda[j][alpha];
}
}
diffsq = New_diffsq;
iter++;
return;
}
*/
/*
void XXZ_gpd_Bethe_State::Iterate_BAE_Newton ()
{
// does one step of a Newton method on the rapidities...
Vect_DP RHSBAE (0.0, base.Nraptot); // contains RHS of BAEs
Vect_CX dlambda (0.0, base.Nraptot); // contains delta lambda computed from Newton's method
SQMat_CX Gaudin (0.0, base.Nraptot);
Vect_INT indx (base.Nraptot);
DP sumtheta = 0.0;
DP arg = 0.0;
DP fn_arg = 0.0;
DP olddiffsq = diffsq;
// Compute the RHS of the BAEs:
int index = 0;
(*this).Compute_tanlambda();
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
sumtheta = 0.0;
for (int k = 0; k < chain.Nstrings; ++k) {
for (int beta = 0; beta < base[k]; ++beta)
if ((chain.Str_L[j] == 1) && (chain.Str_L[k] == 1)) {
sumtheta += atan ((tanlambda[j][alpha] - tanlambda[k][beta])/((1.0 + tanlambda[j][alpha] * tanlambda[k][beta])
* chain.ta_n_anis_over_2[2]))
+ PI * floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
}
else sumtheta += 0.5 * Theta_XXZ_gpd((tanlambda[j][alpha] - tanlambda[k][beta])/(1.0 + tanlambda[j][alpha] * tanlambda[k][beta]),
chain.Str_L[j], chain.Str_L[k], chain.ta_n_anis_over_2)
+ PI * (2.0 * JSC::min(chain.Str_L[j], chain.Str_L[k]) - ((j == k) ? 1.0 : 0))
* floor(0.5 + (lambda[j][alpha] - lambda[k][beta])/PI);
}
sumtheta *= 2.0;
RHSBAE[index] = chain.Nsites * 2.0 * (atan(tanlambda[j][alpha]/chain.ta_n_anis_over_2[chain.Str_L[j]])
+ PI * floor(0.5 + lambda[j][alpha]/PI))
// )
- sumtheta - PI*Ix2[j][alpha];
index++;
}
}
(*this).Build_Reduced_Gaudin_Matrix (Gaudin);
for (int i = 0; i < base.Nraptot; ++i) dlambda[i] = - RHSBAE[i];
DP d;
ludcmp_CX (Gaudin, indx, d);
lubksb_CX (Gaudin, indx, dlambda);
diffsq = 0.0;
// for (int i = 0; i < base.Nraptot; ++i) diffsq += norm(dlambda[i]);
int ctr = 0;
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
diffsq += norm(tan(lambda[j][alpha] + dlambda[ctr]) - tanlambda[j][alpha]);
// cout << "lambda = " << lambda[j][alpha] << "\tdlambda = " << dlambda[ctr] << endl;
ctr++;
}
}
// if we've converged, calculate the norm here, since the work has been done...
if (diffsq < chain.prec) {
lnnorm = 0.0;
for (int j = 0; j < base.Nraptot; j++) lnnorm += log(abs(Gaudin[j][j]));
}
index = 0;
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
lambda[j][alpha] += real(dlambda[index]);
index++;
}
}
iter_Newton++;
// cout << "iter_N = " << iter_Newton << "\t" << diffsq << endl;
return;
}
*/
bool XXZ_gpd_Bethe_State::Check_Rapidities()
{
bool nonan = true;
for (int j = 0; j < chain.Nstrings; ++j)
for (int alpha = 0; alpha < base[j]; ++alpha) if (nonan) nonan = ((!is_nan(lambda[j][alpha]))
//&& (lambda[j][alpha] > -0.5*PI*chain.Str_L[j])
//&& (lambda[j][alpha] <= 0.5*PI*chain.Str_L[j])
);
bool all_within_pi_interval = true;
DP min_lambda = 10.0;
DP max_lambda = -10.0;
for (int j = 0; j < chain.Nstrings; ++j)
for (int alpha = 0; alpha < base[j]; ++alpha) {
if (lambda[j][alpha] < min_lambda) min_lambda = lambda[j][alpha];
if (lambda[j][alpha] > max_lambda) max_lambda = lambda[j][alpha];
}
//all_within_pi_interval = (fabs(max_lambda - min_lambda) < PI - 0.0 * 1.0e-12);
// Use pi interval, but centered on 0:
all_within_pi_interval = max_lambda < 0.5* PI + 1.0e-8 && min_lambda > -0.5* PI - 1.0e-8;
//nonan *= all_within_pi_interval;
return nonan;
}
DP XXZ_gpd_Bethe_State::String_delta ()
{
// Computes the sum of absolute value of \delta^{a, a+1} in string hypothesis, for a given bethe eigenstate
DP delta = 0.0;
int occupied_strings = 0;
for (int i = 0; i < (*this).chain.Nstrings; ++i) if ((*this).chain.Str_L[i] > 1) occupied_strings += (*this).base.Nrap[i];
//if ((*this).conv == 0) delta = 1.0;
if (occupied_strings == 0) delta = 0.0;
else {
Vect_DP ln_deltadiff(0.0, 1000); // contains ln |delta^{a, a+1}|
Vect_DP deltadiff(0.0, 1000); // contains |delta^{a, a+1}|
complex<DP> log_BAE_reg = 0.0;
for (int j = 0; j < (*this).chain.Nstrings; ++j) {
for (int alpha = 0; alpha < (*this).base[j]; ++alpha) {
ln_deltadiff = 0.0;
for (int a = 1; a <= (*this).chain.Str_L[j]; ++a) {
if ((*this).chain.Str_L[j] > 1) { // else the BAE are already 1
log_BAE_reg = DP((*this).chain.Nsites) * log(sin((*this).lambda[j][alpha] + 0.5 * II * (*this).chain.anis
* ((*this).chain.Str_L[j] + 1.0 - 2.0 * a + 1.0))
/sin((*this).lambda[j][alpha] + 0.5 * II * (*this).chain.anis * ((*this).chain.Str_L[j] + 1.0 - 2.0 * a - 1.0)));
for (int k = 0; k < (*this).chain.Nstrings; ++k)
for (int beta = 0; beta < (*this).base[k]; ++beta)
for (int b = 1; b <= (*this).chain.Str_L[k]; ++b) {
if ((j != k) || (alpha != beta) || (a != b - 1))
log_BAE_reg += log(sin(((*this).lambda[j][alpha] + 0.5 * II * (*this).chain.anis * ((*this).chain.Str_L[j] + 1.0 - 2.0 * a ))
- ((*this).lambda[k][beta] + 0.5 * II * (*this).chain.anis * ((*this).chain.Str_L[k] + 1.0 - 2.0 * b ))
- II * (*this).chain.anis));
if ((j != k) || (alpha != beta) || (a != b + 1))
log_BAE_reg -= log(sin(((*this).lambda[j][alpha] + 0.5 * II * (*this).chain.anis * ((*this).chain.Str_L[j] + 1.0 - 2.0 * a ))
- ((*this).lambda[k][beta] + 0.5 * II * (*this).chain.anis * ((*this).chain.Str_L[k] + 1.0 - 2.0 * b ))
+ II * (*this).chain.anis));
}
// The regular LHS of BAE is now defined. Now sum up the deltas...
if (a == 1) ln_deltadiff[0] = -real(log_BAE_reg);
else if (a < (*this).chain.Str_L[j]) ln_deltadiff[a - 1] = ln_deltadiff[a-2] - real(log_BAE_reg);
else if (a == (*this).chain.Str_L[j]) ln_deltadiff[a-1] = real(log_BAE_reg);
} // if ((*this).chain.Str_L[j] > 1)
} // for (int a = 1; ...
for (int a = 0; a < (*this).chain.Str_L[j]; ++a) {
deltadiff[a] = ln_deltadiff[a] != 0.0 ? exp(ln_deltadiff[a]) : 0.0;
delta += fabs(deltadiff[a]);
}
} // alpha sum
} // j sum
if (is_nan(delta)) delta = 1.0; // sentinel
} // else
return delta;
}
void XXZ_gpd_Bethe_State::Compute_Energy ()
{
DP sum = 0.0;
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
sum += sinh(chain.Str_L[j] * chain.anis) / (cos(2.0 * lambda[j][alpha]) - cosh(chain.Str_L[j] * chain.anis));
}
}
sum *= chain.J * sinh(chain.anis);
E = sum;
return;
}
/*
void XXZ_gpd_Bethe_State::Compute_Momentum ()
{
int sum_Ix2 = 0;
DP sum_M = 0.0;
for (int j = 0; j < chain.Nstrings; ++j) {
sum_M += base[j];
for (int alpha = 0; alpha < base[j]; ++alpha) {
sum_Ix2 += Ix2[j][alpha];
}
}
iK = (chain.Nsites/2) * int(sum_M) - (sum_Ix2/2);
while (iK >= chain.Nsites) iK -= chain.Nsites;
while (iK < 0) iK += chain.Nsites;
K = PI * sum_M - PI * sum_Ix2/chain.Nsites;
while (K >= 2.0*PI) K -= 2.0*PI;
while (K < 0.0) K += 2.0*PI;
return;
}
*/
void XXZ_gpd_Bethe_State::Build_Reduced_Gaudin_Matrix (SQMat<complex<DP> >& Gaudin_Red)
{
if (Gaudin_Red.size() != base.Nraptot) JSCerror("Passing matrix of wrong size in Build_Reduced_Gaudin_Matrix.");
int index_jalpha;
int index_kbeta;
DP sum_hbar_XXZ = 0.0;
DP sinhetasq = pow(sinh(chain.anis), 2.0);
(*this).Compute_sinlambda();
(*this).Compute_coslambda();
index_jalpha = 0;
for (int j = 0; j < chain.Nstrings; ++j) {
for (int alpha = 0; alpha < base[j]; ++alpha) {
index_kbeta = 0;
for (int k = 0; k < chain.Nstrings; ++k) {
for (int beta = 0; beta < base[k]; ++beta) {
if ((j == k) && (alpha == beta)) {
sum_hbar_XXZ = 0.0;
for (int kp = 0; kp < chain.Nstrings; ++kp) {
for (int betap = 0; betap < base[kp]; ++betap) {
if (!((j == kp) && (alpha == betap)))
sum_hbar_XXZ
+= ddlambda_Theta_XXZ_gpd (lambda[j][alpha] - lambda[kp][betap], chain.Str_L[j], chain.Str_L[kp],
chain.si_n_anis_over_2);
}
}
Gaudin_Red[index_jalpha][index_kbeta]
= complex<DP> ( chain.Nsites * hbar_XXZ_gpd (lambda[j][alpha], chain.Str_L[j], chain.si_n_anis_over_2) - sum_hbar_XXZ);
}
else {
if ((chain.Str_L[j] == 1) && (chain.Str_L[k] == 1))
Gaudin_Red[index_jalpha][index_kbeta] =
complex<DP> (chain.si_n_anis_over_2[4]/(pow(sinlambda[j][alpha] * coslambda[k][beta]
- coslambda[j][alpha] * sinlambda[k][beta], 2.0) + sinhetasq));
else
Gaudin_Red[index_jalpha][index_kbeta] = complex<DP> (ddlambda_Theta_XXZ_gpd (lambda[j][alpha] - lambda[k][beta], chain.Str_L[j],
chain.Str_L[k], chain.si_n_anis_over_2));
}
index_kbeta++;
}
}
index_jalpha++;
}
}
return;
}
// ****************************************************************************************************
// non-member functions
inline DP fbar_XXZ_gpd (DP tanlambda, DP tanhnetaover2)
{
return (2.0 * atan(tanlambda/tanhnetaover2));
}
DP Theta_XXZ_gpd (DP tanlambda, int nj, int nk, DP* tanhnetaover2)
{
DP result;
if ((nj == 1) && (nk == 1)) result = fbar_XXZ_gpd(tanlambda, tanhnetaover2[2]);
else {
result = (nj == nk) ? 0.0 : fbar_XXZ_gpd(tanlambda, tanhnetaover2[fabs(nj - nk)]);
for (int a = 1; a < JSC::min(nj, nk); ++a) result += 2.0 * fbar_XXZ_gpd(tanlambda, tanhnetaover2[fabs(nj - nk) + 2*a]);
result += fbar_XXZ_gpd(tanlambda, tanhnetaover2[nj + nk]);
}
return (result);
}
DP hbar_XXZ_gpd (DP lambda, int n, DP* si_n_anis_over_2)
{
return (si_n_anis_over_2[2*n]/(pow(sin(lambda), 2.0) + pow(si_n_anis_over_2[n], 2.0)));
}
DP ddlambda_Theta_XXZ_gpd (DP lambda, int nj, int nk, DP* si_n_anis_over_2)
{
DP result = (nj == nk) ? 0.0 : hbar_XXZ_gpd(lambda, fabs(nj - nk), si_n_anis_over_2);
for (int a = 1; a < JSC::min(nj, nk); ++a) result += 2.0 * hbar_XXZ_gpd(lambda, fabs(nj - nk) + 2*a, si_n_anis_over_2);
result += hbar_XXZ_gpd(lambda, nj + nk, si_n_anis_over_2);
return (result);
}
XXZ_gpd_Bethe_State Add_Particle_at_Center (const XXZ_gpd_Bethe_State& RefState)
{
if (2*RefState.base.Mdown == RefState.chain.Nsites)
JSCerror("Trying to add a down spin to a zero-magnetized chain in Add_Particle_at_Center.");
Vect<int> newM = RefState.base.Nrap;
newM[0] = newM[0] + 1;
Heis_Base newBase (RefState.chain, newM);
XXZ_gpd_Bethe_State ReturnState (RefState.chain, newBase);
for (int il = 1; il < RefState.chain.Nstrings; ++il)
for (int alpha = 0; alpha < RefState.base.Nrap[il]; ++alpha)
ReturnState.Ix2[il][alpha] = RefState.Ix2[il][alpha];
// Add a quantum number in middle (explicitly: to right of index M[0]/2)
// and shift quantum numbers by half-integer away from added one:
ReturnState.Ix2[0][RefState.base.Nrap[0]/2] = RefState.Ix2[0][RefState.base.Nrap[0]/2] - 1;
for (int i = 0; i < RefState.base.Nrap[0] + 1; ++i)
ReturnState.Ix2[0][i + (i >= RefState.base.Nrap[0]/2)] = RefState.Ix2[0][i] - 1 + 2*(i >= RefState.base.Nrap[0]/2);
return(ReturnState);
}
XXZ_gpd_Bethe_State Remove_Particle_at_Center (const XXZ_gpd_Bethe_State& RefState)
{
if (RefState.base.Nrap[0] == 0)
JSCerror("Trying to remove a down spin in an empty Nrap[0] state.");
Vect<int> newM = RefState.base.Nrap;
newM[0] = newM[0] - 1;
Heis_Base newBase (RefState.chain, newM);
XXZ_gpd_Bethe_State ReturnState (RefState.chain, newBase);
for (int il = 1; il < RefState.chain.Nstrings; ++il)
for (int alpha = 0; alpha < RefState.base.Nrap[il]; ++alpha)
ReturnState.Ix2[il][alpha] = RefState.Ix2[il][alpha];
// Remove midmost and shift quantum numbers by half-integer towards removed one:
for (int i = 0; i < RefState.base.Nrap[0]-1; ++i)
ReturnState.Ix2[0][i] = RefState.Ix2[0][i + (i >= RefState.base.Nrap[0]/2)] + 1 - 2*(i >= RefState.base.Nrap[0]/2);
return(ReturnState);
}
} // namespace JSC