ABACUS/src/YOUNG/Young_Tableau.cc

871 lines
28 KiB
C++

/**********************************************************
This software is part of J.-S. Caux's ABACUS library.
Copyright (c) J.-S. Caux.
-----------------------------------------------------------
File: Young_Tableau.cc
Purpose: Defines Young_Tableau class and procedures.
***********************************************************/
#include "ABACUS.h"
using namespace std;
namespace ABACUS {
Young_Tableau::Young_Tableau () : Nrows(0), Ncols(0), Row_L(0), Col_L(0), id(0LL), maxid(0LL),
map(new long long int [1]), map_computed(false),
idnr_reached(0LL), nboxes_reached(-1),
dimchoose(0), choose_table(0LL)
{}
Young_Tableau::Young_Tableau (int Nr, int Nc)
: Nrows(Nr), Ncols(Nc), Row_L(new int[Nrows]), Col_L(new int[Ncols]), id(0LL),
maxid(choose_lli(Nr + Nc, Nc) - 1LL),
map(new long long int [YOUNG_TABLEAU_ID_OPTION == 2 ? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1]),
map_computed(false), idnr_reached(0LL), nboxes_reached(-1),
dimchoose (ABACUS::min(Nr, Nc) + 1),
choose_table(new long long int[(Nr + Nc + 1) * dimchoose])
{
// Constructs empty tableau of appropriate size
for (int i = 0; i < Nrows; ++i) Row_L[i] = 0;
for (int i = 0; i < Ncols; ++i) Col_L[i] = 0;
// Construct the choose_table
for (int cti = 0; cti < Nr + Nc + 1; ++cti)
for (int ctj = 0; ctj < dimchoose; ++ctj) {
if (cti >= ctj) choose_table[dimchoose * cti + ctj] = choose_lli(cti, ctj);
else choose_table[dimchoose * cti + ctj] = 0LL;
}
// Fill map with zeros:
for (int i = 0; i < (YOUNG_TABLEAU_ID_OPTION == 2 ? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1); ++i)
map[i] = 0LL;
}
Young_Tableau::Young_Tableau (const Young_Tableau& RefTableau) // copy constructor
: Nrows(RefTableau.Nrows), Ncols(RefTableau.Ncols), Row_L(new int[RefTableau.Nrows]), Col_L(new int[RefTableau.Ncols]),
id(RefTableau.id), maxid(RefTableau.maxid),
map(new long long int [YOUNG_TABLEAU_ID_OPTION == 2 ? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1]),
map_computed(RefTableau.map_computed), idnr_reached(RefTableau.idnr_reached),
nboxes_reached(RefTableau.nboxes_reached),
dimchoose (RefTableau.dimchoose),
choose_table(new long long int[(Nrows + Ncols + 1) * dimchoose])
{
for (int i = 0; i < Nrows; ++i) Row_L[i] = RefTableau.Row_L[i];
for (int i = 0; i < Ncols; ++i) Col_L[i] = RefTableau.Col_L[i];
// Construct the choose_table
for (int cti = 0; cti < Nrows + Ncols + 1; ++cti)
for (int ctj = 0; ctj < dimchoose; ++ctj) {
if (cti >= ctj) choose_table[dimchoose * cti + ctj] = choose_lli(cti, ctj);
else choose_table[dimchoose * cti + ctj] = 0LL;
}
// The map:
for (int i = 0; i < (YOUNG_TABLEAU_ID_OPTION == 2 ? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1); ++i)
map[i] = RefTableau.map[i];
}
Young_Tableau::Young_Tableau (int Nr, int Nc, const Young_Tableau& RefTableau)
: Nrows(Nr), Ncols(Nc), Row_L(new int[Nrows]), Col_L(new int[Ncols]), id(0LL),
maxid(choose_lli(Nr + Nc, Nc) - 1LL),
map(new long long int[YOUNG_TABLEAU_ID_OPTION == 2 ? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1]),
map_computed(false), idnr_reached(0LL), nboxes_reached(-1),
dimchoose (ABACUS::min(Nr, Nc) + 1),
choose_table(new long long int[(Nr + Nc + 1) * dimchoose])
{
// Constructs empty tableau of appropriate size
for (int i = 0; i < Nrows; ++i) Row_L[i] = 0;
for (int i = 0; i < Ncols; ++i) Col_L[i] = 0;
// Construct the choose_table
// Copy entries from reference table
for (int cti = 0; cti < ABACUS::min(Nr + Nc + 1, RefTableau.Nrows + RefTableau.Ncols + 1); ++cti)
for (int ctj = 0; ctj < ABACUS::min(dimchoose, RefTableau.dimchoose); ++ctj)
choose_table[dimchoose * cti + ctj] = cti >= ctj ? RefTableau.choose_table[RefTableau.dimchoose * cti + ctj] : 0LL;
// add missing parts if there are any
int refdim1 = RefTableau.Nrows + RefTableau.Ncols + 1;
for (int cti = 0; cti < Nr + Nc + 1; ++cti)
for (int ctj = 0; ctj < dimchoose; ++ctj)
if (cti >= refdim1 || ctj >= RefTableau.dimchoose)
choose_table[dimchoose * cti + ctj] = cti >= ctj ? choose_lli(cti, ctj) : 0LL;
// The map:
for (int i = 0; i < (YOUNG_TABLEAU_ID_OPTION == 2 ? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1); ++i)
map[i] = 0LL;
}
Young_Tableau& Young_Tableau::operator= (const Young_Tableau& RefTableau)
{
if (this != &RefTableau) {
Nrows = RefTableau.Nrows;
Ncols = RefTableau.Ncols;
if (Row_L != 0) delete[] Row_L;
Row_L = new int[Nrows];
for (int i = 0; i < Nrows; ++i) Row_L[i] = RefTableau.Row_L[i];
if (Col_L != 0) delete[] Col_L;
Col_L = new int[Ncols];
for (int i = 0; i < Ncols; ++i) Col_L[i] = RefTableau.Col_L[i];
id = RefTableau.id;
maxid = RefTableau.maxid;
dimchoose = RefTableau.dimchoose;
if (choose_table != 0LL) delete[] choose_table;
choose_table = new long long int[(Nrows + Ncols + 1) * dimchoose];
for (int cti = 0; cti < Nrows + Ncols + 1; ++cti)
for (int ctj = 0; ctj < dimchoose; ++ctj)
{
choose_table[dimchoose * cti + ctj] = RefTableau.choose_table[dimchoose * cti + ctj];
}
if (map != 0LL) delete[] map;
map = new long long int[YOUNG_TABLEAU_ID_OPTION == 2 ? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1];
for (long long int i = 0; i < (YOUNG_TABLEAU_ID_OPTION == 2
? ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT) : 1); ++i)
map[i] = RefTableau.map[i];
map_computed = RefTableau.map_computed;
idnr_reached = RefTableau.idnr_reached;
nboxes_reached = RefTableau.nboxes_reached;
}
if (choose_table[0] != 1LL) { // compute the table
for (int cti = 0; cti < Nrows + Ncols + 1; ++cti)
for (int ctj = 0; ctj < dimchoose; ++ctj)
{
choose_table[dimchoose * cti + ctj] = cti >= ctj ? choose_lli(cti, ctj) : 0LL;
}
}
return(*this);
}
Young_Tableau::~Young_Tableau () // destructor
{
if (Row_L != 0) delete[] Row_L;
if (Col_L != 0) delete[] Col_L;
if (choose_table != 0) delete[] choose_table;
if (map != 0LL) delete[] map;
}
//*********************************************************************************
// Member functions
Young_Tableau& Young_Tableau::Compute_Map (long long int idnr_to_reach)
{
if (idnr_to_reach < 0LL) ABACUSerror("negative id requested in Compute_Map");
long long int idnr_to_reach_here = ABACUS::min(idnr_to_reach, ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT));
if (!map_computed) {
nboxes_reached = -1;
idnr_reached = 0LL;
map_computed = true;
}
if (idnr_to_reach_here >= idnr_reached) {
if (maxid == 0) map[0] = 0LL;
else {
// Keep a copy of the Tableau, in case...
int* old_Row_L = new int[Nrows];
for (int i = 0; i < Nrows; ++i) old_Row_L[i] = Row_L[i];
while (idnr_reached <= ABACUS::min(idnr_to_reach, TABLEAU_ID_UPPER_LIMIT)) {
nboxes_reached++;
Distribute_boxes (nboxes_reached, 0);
}
// Reset Row_L
for (int i = 0; i < Nrows; ++i) Row_L[i] = old_Row_L[i];
delete[] old_Row_L;
Set_Col_L_given_Row_L();
}
}
return(*this);
}
Young_Tableau& Young_Tableau::Distribute_boxes (int nboxes_to_dist, int level)
{
if (level > Nrows) {
// Failed to achieve proper distribution. Do nothing.
}
else if (nboxes_to_dist > Ncols * (Nrows - level)) {
cout << Nrows << "\t" << Ncols << "\t" << level << "\t" << nboxes_to_dist << "\t"
<< idnr_reached << "\t" << nboxes_reached << endl;
ABACUSerror("nboxes_to_dist too high");
}
else if (nboxes_to_dist == 0) {
for (int j = level; j < Nrows; ++j) Row_L[j] = 0;
Compute_id (0);
if (idnr_reached < TABLEAU_ID_UPPER_LIMIT) map[idnr_reached] = id;
idnr_reached++;
}
else {
int maxlength = level == 0 ? Ncols : Row_L[level - 1];
for (int nboxes_of_level = (nboxes_to_dist - 1)/(Nrows - level) + 1;
nboxes_of_level <= ABACUS::min(nboxes_to_dist, maxlength); ++nboxes_of_level) {
Row_L[level] = nboxes_of_level;
Distribute_boxes (nboxes_to_dist - nboxes_of_level, level + 1);
}
}
return(*this);
}
Young_Tableau& Young_Tableau::Compute_id()
{
Compute_id(0);
return(*this);
}
long long int Young_Tableau::Compute_Descendent_id (int option, Vect_INT& Desc_Row_L, int Nrows_Desc, int Ncols_Desc,
const Young_Tableau& RefTableau)
{
long long int answer = 0;
if (option == 0) {
if (Nrows_Desc == 0) answer = 0;
else if (Desc_Row_L[0] == 0) answer = 0;
else if (Nrows_Desc == 1) answer = Desc_Row_L[0];
else {
int highest_occupied_row = Nrows_Desc - 1;
while (Desc_Row_L[highest_occupied_row] == 0) highest_occupied_row--; // index of highest occupied row;
for (int j = 0; j < Desc_Row_L[highest_occupied_row]; ++j)
answer += RefTableau.choose_table[RefTableau.dimchoose * (highest_occupied_row + Ncols_Desc - j)
+ ABACUS::min(highest_occupied_row, Ncols_Desc - j)];
Vect_INT Desc_Desc_Row_L(highest_occupied_row);
for (int i = 0; i < highest_occupied_row; ++i)
Desc_Desc_Row_L[i] = Desc_Row_L[i] - Desc_Row_L[highest_occupied_row];
answer += Compute_Descendent_id (0, Desc_Desc_Row_L, highest_occupied_row,
Ncols_Desc - Desc_Row_L[highest_occupied_row], RefTableau);
}
}
else if (option == 1) {
if (Nrows_Desc == 0 || Ncols_Desc == 0) answer = 0;
else if (Desc_Row_L[0] == 0) answer = 0;
else if ((Nrows_Desc == 1) || (Desc_Row_L[1] == 0)) answer = Desc_Row_L[0];
else {
int ndiag = 0;
while (Desc_Row_L[ndiag] > ndiag) ndiag++;
if (ndiag == 1) {
int Desc_Col_L_0 = Nrows_Desc;
while (Desc_Row_L[Desc_Col_L_0 - 1] == 0) Desc_Col_L_0--;
answer = (Desc_Col_L_0 - 1) * Ncols_Desc + Desc_Row_L[0];
}
else {
for (int j = 0; j < ndiag; ++j) answer += RefTableau.choose_table[RefTableau.dimchoose * Nrows_Desc + j]
* RefTableau.choose_table[RefTableau.dimchoose * Ncols_Desc + j];
Vect_INT Desc1_Row_L(ndiag);
for (int i = 0; i < ndiag; ++i) Desc1_Row_L[i] = Desc_Row_L[i] - ndiag;
Vect_INT Desc2_Row_L(ndiag);
for (int i = 0; i < ndiag; ++i) {
Desc2_Row_L[i] = 0;
while(Desc_Row_L[Desc2_Row_L[i] ] > i && Desc2_Row_L[i] < Ncols_Desc) Desc2_Row_L[i]++;
Desc2_Row_L[i] -= ndiag;
}
long long int Desc1_id = Compute_Descendent_id (1, Desc1_Row_L, ndiag, Ncols_Desc - ndiag, RefTableau);
long long int Desc2_id = Compute_Descendent_id (1, Desc2_Row_L, ndiag, Nrows_Desc - ndiag, RefTableau);
answer += Desc2_id * RefTableau.choose_table[RefTableau.dimchoose * Ncols_Desc + ndiag] + Desc1_id;
}
}
} // else if (option == 1)
else ABACUSerror("Wrong option in Young_Tableau::Compute_Descendent_id");
return(answer);
}
Young_Tableau& Young_Tableau::Compute_id (int option) // computes the id according to rule option
{
long long int idnr = 0;
if (option == 0) {
if (Nrows == 0) idnr = 0;
else if (Row_L[0] == 0) idnr = 0;
else if (Nrows == 1) idnr = Row_L[0];
else {
int highest_occupied_row = Nrows - 1;
while (Row_L[highest_occupied_row] == 0) highest_occupied_row--; // index of highest occupied row;
for (int j = 0; j < Row_L[highest_occupied_row]; ++j)
idnr += choose_table[dimchoose * (highest_occupied_row + Ncols - j)
+ ABACUS::min(highest_occupied_row, Ncols - j)];
Vect_INT Desc_Row_L(highest_occupied_row);
for (int i = 0; i < highest_occupied_row; ++i) Desc_Row_L[i] = Row_L[i] - Row_L[highest_occupied_row];
idnr += Compute_Descendent_id (0, Desc_Row_L, highest_occupied_row, Ncols - Row_L[highest_occupied_row], (*this));
}
}
else if (option == 1) {
if (Nrows == 0 || Ncols == 0) idnr = 0;
else if (Row_L[0] == 0) idnr = 0;
else if ((Nrows == 1) || (Row_L[1] == 0)) idnr = Row_L[0];
else {
int ndiag = 0;
while (Row_L[ndiag] > ndiag) ndiag++;
if (ndiag == 1) {
idnr = (Col_L[0] - 1) * Ncols + Row_L[0];
}
else {
for (int j = 0; j < ndiag; ++j) idnr += choose_table[dimchoose * Nrows + j] * choose_table[dimchoose * Ncols + j];
Vect_INT Desc1_Row_L(ndiag);
for (int i = 0; i < ndiag; ++i) Desc1_Row_L[i] = Row_L[i] - ndiag;
Vect_INT Desc2_Row_L(ndiag);
for (int i = 0; i < ndiag; ++i) Desc2_Row_L[i] = Col_L[i] - ndiag;
long long int Desc1_id = Compute_Descendent_id (1, Desc1_Row_L, ndiag, Ncols - ndiag, (*this));
long long int Desc2_id = Compute_Descendent_id (1, Desc2_Row_L, ndiag, Nrows - ndiag, (*this));
idnr += Desc2_id * choose_table[dimchoose * Ncols + ndiag] + Desc1_id;
}
}
} // else if (option == 1)
else if (option == 2) {
// The order here is given by first listing all tableaux with 0 boxes, then 1, 2, ... by using Map
Compute_id (0); // sets the id according to rule 0
Compute_Map (idnr); // make sure the state map is computed
while (map[idnr] != id && idnr < ABACUS::min(maxid + 1LL, TABLEAU_ID_UPPER_LIMIT))
idnr++; // match with inverse map to get the idnr according to rule 2
}
else ABACUSerror("Wrong option for Tableau ids");
id = idnr;
return (*this);
}
Young_Tableau& Young_Tableau::Set_to_id (long long int idnr)
{
Set_to_id (idnr, 0);
return(*this);
}
Young_Tableau& Young_Tableau::Set_to_id (long long int idnr, int option)
// sets the tableau to the one corresponding to idnr
{
if (option == 0) {
if ((idnr < 0) || ((maxid < idnr) && (Nrows*Ncols != 0))) {
cout << "Nrows = " << Nrows << "\tNcols = " << Ncols
<< "\tmaxid = " << maxid << "\trequested id = " << idnr << endl;
ABACUSerror("Wrong idnr in Set_to_id for Young Tableau.");
}
id = idnr;
int NColumnseff = Ncols;
long long int ideff = idnr;
for (int i = 0; i < Ncols; ++i) {
if (ideff == 0) Col_L[i] = 0;
else {
int nbar = 0;
if (ideff <= 0) nbar = 0;
while (ideff >= choose_table[dimchoose * (NColumnseff + nbar) + ABACUS::min(NColumnseff, nbar)]) nbar++;
Col_L[i] = nbar;
ideff -= choose_table[dimchoose * (NColumnseff + Col_L[i] - 1) + ABACUS::min(NColumnseff, Col_L[i] - 1)];
NColumnseff--;
}
} // for (int i = 0; i < Ncols...
(*this).Set_Row_L_given_Col_L();
} // if (option == 0)
else if (option == 1) {
if ((idnr < 0LL) || ((maxid < idnr) && (Nrows*Ncols != 0)))
ABACUSerror("Wrong idnr in Set_to_id for Young Tableau.");
if (Nrows*Ncols == 0 && idnr != 0LL) ABACUSerror("Trying nonzero id on empty Tableau.");
id = idnr;
if (idnr == 0LL) {
for (int i = 0; i < Nrows; ++i) Row_L[i] = 0;
for (int i = 0; i < Ncols; ++i) Col_L[i] = 0;
}
else {
int ndiag = 0;
int sum = 0;
while ((sum < idnr) && (ndiag < Nrows) && (ndiag < Ncols)) {
ndiag++;
sum += choose_table[dimchoose * Nrows + ndiag] * choose_table[dimchoose * Ncols + ndiag];
}
long long int residual_id = idnr - 1 - sum + choose_table[dimchoose * Nrows + ndiag]
* choose_table[dimchoose * Ncols + ndiag];
if (ndiag == 0 && idnr != 0LL) ABACUSerror("Zero ndiag for nonzero idnr in Tableau.");
else if (ndiag == 1) {
if (Nrows >= Ncols) {
Col_L[0] = (idnr - 1)/Ncols + 1;
for (int i = 1; i < idnr - Ncols * (Col_L[0] - 1); ++i) Col_L[i] = 1;
for (int i = idnr - Ncols * (Col_L[0] - 1); i < Ncols; ++i) Col_L[i] = 0;
(*this).Set_Row_L_given_Col_L();
}
else if (Nrows < Ncols) {
Row_L[0] = (idnr - 1)/Nrows + 1;
for (int i = 1; i < idnr - Nrows * (Row_L[0] - 1); ++i) Row_L[i] = 1;
for (int i = idnr - Nrows * (Row_L[0] - 1); i < Nrows; ++i) Row_L[i] = 0;
(*this).Set_Col_L_given_Row_L();
}
}
else {
Young_Tableau Residual1(ndiag, Ncols - ndiag);
Young_Tableau Residual2(ndiag, Nrows - ndiag);
Residual2.Set_to_id(residual_id/(Residual1.maxid + 1LL));
residual_id -= (Residual1.maxid + 1LL) * Residual2.id;
Residual1.Set_to_id(residual_id);
for (int i = 0; i < ndiag; ++i) Row_L[i] = ndiag + Residual1.Row_L[i];
for (int i = 0; i < Residual2.Ncols; ++i) Row_L[ndiag + i] = Residual2.Col_L[i];
(*this).Set_Col_L_given_Row_L();
}
}
} // else if (option == 1)
else if (option == 2) {
Compute_Map (idnr); // make sure the state map is computed
(*this).Set_to_id (map[idnr], 0);
id = idnr;
}
else ABACUSerror("Wrong option for Tableau id in Set_to_id");
return(*this);
}
Young_Tableau& Young_Tableau::Set_Row_L (Vect_INT& Row_Lengths) // set row lengths to elements of given vector
{
if (Row_Lengths.size() != Nrows)
ABACUSerror("Vector of incompatible dimension used to initialize Young Tableau.");
for (int i = 0; i < Row_Lengths.size() - 1; ++i)
if (Row_Lengths[i] < Row_Lengths[i+1])
ABACUSerror("Vector is not a proper Young tableau.");
for (int i = 0; i < Nrows; ++i) Row_L[i] = Row_Lengths[i];
(*this).Set_Col_L_given_Row_L();
return(*this);
}
Young_Tableau& Young_Tableau::Set_Col_L_given_Row_L () // sets the Col_L array self-consistently
{
int col_l = 0;
for (int i = 0; i < Ncols; ++i) {
col_l = 0;
while ((Row_L[col_l] > i) && (col_l < Nrows)) ++col_l;
Col_L[i] = col_l;
}
return *this;
}
Young_Tableau& Young_Tableau::Set_Row_L_given_Col_L () // sets the Col_L array self-consistently
{
int row_l = 0;
for (int i = 0; i < Nrows; ++i) {
row_l = 0;
while ((Col_L[row_l] > i) && (row_l < Ncols)) ++row_l;
Row_L[i] = row_l;
}
return *this;
}
Young_Tableau& Young_Tableau::Print () // couts the tableau
{
cout << endl;
for (int i = 0; i < Nrows; ++i) {
if (Row_L[i] > 0) {
for (int j = 0; j < Row_L[i]; ++j) cout << "X";
cout << endl;
}
}
cout << endl;
return(*this);
}
std::ostream& operator<< (std::ostream& s, const Young_Tableau& tableau)
{
s << endl;
for (int i = 0; i < tableau.Nrows; ++i) {
if (tableau.Row_L[i] > 0) {
for (int j = 0; j < tableau.Row_L[i]; ++j) s << "X";
s << endl;
}
}
s << endl;
return(s);
}
bool Young_Tableau::Lower_Row (int i)
{
// Removes a box from row i. Returns success boolean.
// Recomputes id.
if (id == 0LL || Nrows == 0 || Ncols == 0) return(false); // Tableau is empty
if (i < 0 || i >= Nrows) ABACUSerror("Trying to Lower_Row of out of bounds index.");
if (Row_L[i] < 1 || (i < Nrows - 1 && Row_L[i+1] == Row_L[i])) // can't lower Row_L, breaks tableau rules
return(false);
else {
Row_L[i] -= 1;
Set_Col_L_given_Row_L();
Compute_id();
return(true);
}
return(false);
}
bool Young_Tableau::Raise_Row (int i)
{
// Adds a box to row i. Returns success boolean.
// Recomputes id.
if (Nrows == 0 || Ncols == 0) return(false); // Tableau is empty
if (i < 0 || i >= Nrows) ABACUSerror("Trying to Raise_Row of out of bounds index.");
if (Row_L[i] == Ncols || (i > 1 && Row_L[i-1] == Row_L[i])) // can't raise Row_L, breaks tableau rules
return(false);
else {
Row_L[i] += 1;
Set_Col_L_given_Row_L();
Compute_id();
return(true);
}
return(false);
}
bool Young_Tableau::Raise_Lowest_Nonzero_Row()
{
// adds a box to the lowest nonzero length Row, recomputes id, returns true if tableau has changed
if (id == 0LL || Nrows == 0 || Ncols == 0) return(false); // Tableau is empty
// otherwise find the lowest nonzero row:
int iln0r = Nrows - 1;
while (Row_L[iln0r] == 0 && iln0r >= 0) iln0r--;
if (iln0r < 0) ABACUSerror("id wrongly set in Young_Tableau (Raise_Lowest_Nonzero_Row).");
// This should not happen, since if iln0r == -1, id should be 0.
else if (iln0r == 0 && Row_L[0] < Ncols || iln0r > 0 && Row_L[iln0r - 1] > Row_L[iln0r]) {
// there is space for at least one more box !
Row_L[iln0r] += 1;
Set_Col_L_given_Row_L();
Compute_id();
return(true);
}
return(false);
}
bool Young_Tableau::Raise_Next_to_Lowest_Nonzero_Row()
{
// same thing, but for Row under lowest nonzero length one.
// Important: allow raising first row if tableau is empty.
if (Ncols == 0 || Nrows == 0) return(false); // no space !
// Find index of lowest nonzero row: can be -1 if Tableau is empty
int iln0r = Nrows - 1;
while (Row_L[iln0r] == 0 && iln0r >= 0) iln0r--;
if (iln0r == -1 && Row_L[0] < Ncols || iln0r >= 0 && iln0r < Nrows - 1 && Row_L[iln0r] > Row_L[iln0r + 1]) {
// there is space for at least one more box !
Row_L[iln0r + 1] += 1;
Set_Col_L_given_Row_L();
Compute_id();
return(true);
}
return(false);
}
bool Young_Tableau::Move_Box_from_Col_to_Col (int ifrom, int ito)
{
// Moves a box from column ifrom to column ito. Recomputes id.
// If fails: returns false, leaves Tableau unchanged.
if (!(ifrom >= 0 && ifrom < Ncols && ito >= 0 && ito < Ncols))
return(false);
else { // try it out
int* Col_L_check = new int[Ncols];
for (int i = 0; i < Ncols; ++i) Col_L_check[i] = Col_L[i];
Col_L_check[ifrom] -= 1;
Col_L_check[ito] += 1;
if (Col_L_check[ifrom] < 0 || Col_L_check[ito] > Nrows) return(false);
for (int i = 0; i < Ncols - 1; ++i)
if (Col_L_check[i] < Col_L_check[i+1]) return(false);
// If we're here, it works.
Col_L[ifrom] -= 1;
Col_L[ito] += 1;
Set_Row_L_given_Col_L();
Compute_id();
}
return(true);
}
Vect<Young_Tableau> Young_Tableau::Descendents (int fixed_Nboxes)
{
// Produces a vector of Young_Tableau which are the descendents of
// the (*this) object. If the number of boxes is not fixed, there
// are up to 2 decendents. If the number of boxes is fixed, there
// can be many more.
int ndesc = 0;
if (!fixed_Nboxes) {
// Here, we can descend by adding a box to the lowest nonzero row,
// to to the empty one immediately below the lowest nonzero row:
Young_Tableau Tableau_check1 = (*this);
bool check1 = (Tableau_check1.Raise_Lowest_Nonzero_Row());
if (check1) ndesc++;
Young_Tableau Tableau_check2 = (*this);
bool check2 = (Tableau_check2.Raise_Next_to_Lowest_Nonzero_Row());
if (check2) ndesc++;
Vect<Young_Tableau> Tableau_desc(ndesc);
if (check1) Tableau_desc[0] = Tableau_check1;
if (check2) Tableau_desc[check1] = Tableau_check2;
return(Tableau_desc);
}
else if (fixed_Nboxes) {
// Here, we can promote a box from the right to the left,
// provided it is in or to the left of col_bdry, where col_bdry
// is the column index of the leftmost column having a box which
// isn't `shadowed' when superimposing the minimal tableau having this
// number of boxes (all boxes in highest possible rows) with the
// actual tableau.
// Define the original lowest tableau:
Young_Tableau Tableau_zero(*this);
int Nboxes_tot = 0;
for (int i = 0; i < Nrows; ++i) Nboxes_tot += Row_L[i];
int Nboxes_above = 0;
for (int level = 0; level < Nrows; ++level) {
Tableau_zero.Row_L[level] = ABACUS::min(Nboxes_tot - Nboxes_above, Ncols);
Nboxes_above += Tableau_zero.Row_L[level];
}
Tableau_zero.Set_Col_L_given_Row_L();
int level_bdry = 0;
for (int level = 0; level < Nrows; ++level)
if (Row_L[level] < Tableau_zero.Row_L[level]) level_bdry = level;
int right_bdry = ABACUS::min(Row_L[level_bdry] + 1, Ncols - 1);
// We can now displace a box from right to left, starting from
// a column with index <= right_bdry.
int left_bdry = 0;
for (int level = 0; level < Ncols; ++level)
if (Col_L[level] > Tableau_zero.Col_L[level]) left_bdry = level;
// We can put a box into a column with index >= left_bdry
// Now do the descendents:
Vect<Young_Tableau> Tableau_desc_init ((right_bdry - left_bdry) * (right_bdry - left_bdry));
Young_Tableau Tableau_ref = (*this);
Young_Tableau Tableau_check = (*this);
for (int ifrom = right_bdry; ifrom >= left_bdry + 1; --ifrom)
for (int ito = ifrom - 1; ito >= left_bdry; --ito) {
if (Tableau_check.Move_Box_from_Col_to_Col (ifrom, ito)) {
Tableau_desc_init[ndesc++] = Tableau_check;
Tableau_check = Tableau_ref;
}
}
Vect<Young_Tableau> Tableau_desc (ndesc);
for (int i = 0; i < ndesc; ++i) Tableau_desc[i] = Tableau_desc_init[i];
return(Tableau_desc);
} // if (fixed_Nboxes)
return(Vect<Young_Tableau> (0));
}
Vect<Young_Tableau> Young_Tableau::Descendents_Boosted_State (int fixed_Nboxes)
{
// Produces a vector of Young_Tableau which are the descendents of
// the (*this) object. If the number of boxes is not fixed, there
// are up to 2 decendents. If the number of boxes is fixed, there
// can be many more.
// IMPORTANT ASSUMPTIONS:
// (*this) is (or is descended from) a boosted state. The only
// descendents considered here are thus those for which the raised
// box is the highest still occupied box of the originally boosted state.
int ndesc = 0;
// Is tableau non-empty ?
if (Nrows <= 0 || Ncols <= 0) return(Vect<Young_Tableau> (0));
int Nboxes = 0;
for (int i = 0; i < Nrows; ++i) Nboxes += Row_L[i];
if (Nboxes == 0) return(Vect<Young_Tableau> (0));
// Look for the level with the highest as yet unraised box of initial
// boosted state:
int level_from = 0;
if (id == maxid) level_from = Nrows - 1; // full tableau treated here
while (Row_L[level_from] == Ncols) level_from++; // necessarily is a row length < Ncols
if (Row_L[level_from] == 0) level_from--; // if row empty or beyond limit, go back one
if (!fixed_Nboxes) {
// The convention here is that we *remove* the highest yet unraised box only
Young_Tableau descendent_attempt = (*this);
if (descendent_attempt.Lower_Row(level_from)) ndesc = 1;
//if (ndesc > 0) {
if (ndesc == 1) {
Vect<Young_Tableau> Tableau_desc(ndesc);
Tableau_desc[0] = descendent_attempt;
return(Tableau_desc);
}
else if (ndesc != 0)
ABACUSerror("There should be either 0 or 1 descendents in Descended_Boosted_State with fixed_iK == true.");
} // if (!fixed_Nboxes)
else if (fixed_Nboxes) {
// We here move the same highest unraised box as for !fixed_Nboxes,
// except that we put higher up in this current tableau:
// We already know from which level we take the box.
// We put it back either on the lowest level, or next:
Young_Tableau Tableau_ref = (*this);
Young_Tableau Tableau_check1 = (*this);
bool check1 = (Tableau_check1.Lower_Row(level_from) && Tableau_check1.Raise_Lowest_Nonzero_Row());
if (check1 && Tableau_check1.Row_L[level_from] == Tableau_ref.Row_L[level_from] - 1)
ndesc++; // to make sure we don't Raise the one we've just removed
Young_Tableau Tableau_check2 = (*this);
bool check2 = (Tableau_check2.Lower_Row(level_from) && Tableau_check2.Raise_Next_to_Lowest_Nonzero_Row());
if (check2 && Tableau_check2.Row_L[level_from] == Tableau_ref.Row_L[level_from] - 1)
ndesc++; // to make sure we don't Raise the one we've just removed
if (ndesc > 0) {
Vect<Young_Tableau> Tableau_desc(ndesc);
if (check1) Tableau_desc[0] = Tableau_check1;
if (check2) Tableau_desc[ndesc - 1] = Tableau_check2; // only up to 2 descendents
return(Tableau_desc);
}
} // if (fixed_Nboxes)
return(Vect<Young_Tableau> (0));
}
int Young_Tableau::Add_Boxes_From_Lowest (int Nboxes)
{
// tries to add Nboxes to Tableau, returns number of boxes added.
if (Ncols == 0 || Nrows == 0) return(0); // can't do anything !
int Nboxes_added = 0;
int previous_Row_L = 0;
for (int working_level = 0; working_level < Nrows; ++working_level) {
previous_Row_L = Row_L[working_level];
Row_L[working_level] = ABACUS::min(previous_Row_L + Nboxes - Nboxes_added, Ncols);
Nboxes_added += Row_L[working_level] - previous_Row_L;
if (Nboxes_added == Nboxes) break;
}
Set_Col_L_given_Row_L();
Compute_id();
return(Nboxes_added);
}
} // namespace ABACUS