#include "TTreeIndex.h"
#include "TTree.h"
#include "TMath.h"
ClassImp(TTreeIndex)
struct IndexSortComparator {
IndexSortComparator(Long64_t *major, Long64_t *minor)
: fValMajor(major), fValMinor(minor)
{}
template<typename Index>
bool operator()(Index i1, Index i2) {
if( *(fValMajor + i1) == *(fValMajor + i2) )
return *(fValMinor + i1) < *(fValMinor + i2);
else
return *(fValMajor + i1) < *(fValMajor + i2);
}
Long64_t *fValMajor, *fValMinor;
};
TTreeIndex::TTreeIndex(): TVirtualIndex()
{
fTree = 0;
fN = 0;
fIndexValues = 0;
fIndexValuesMinor = 0;
fIndex = 0;
fMajorFormula = 0;
fMinorFormula = 0;
fMajorFormulaParent = 0;
fMinorFormulaParent = 0;
}
TTreeIndex::TTreeIndex(const TTree *T, const char *majorname, const char *minorname)
: TVirtualIndex()
{
fTree = (TTree*)T;
fN = 0;
fIndexValues = 0;
fIndexValuesMinor = 0;
fIndex = 0;
fMajorFormula = 0;
fMinorFormula = 0;
fMajorFormulaParent = 0;
fMinorFormulaParent = 0;
fMajorName = majorname;
fMinorName = minorname;
if (!T) return;
fN = T->GetEntries();
if (fN <= 0) {
MakeZombie();
Error("TreeIndex","Cannot build a TreeIndex with a Tree having no entries");
return;
}
GetMajorFormula();
GetMinorFormula();
if (!fMajorFormula || !fMinorFormula) {
MakeZombie();
Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
return;
}
if ((fMajorFormula->GetNdim() != 1) || (fMinorFormula->GetNdim() != 1)) {
MakeZombie();
Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
return;
}
Long64_t *tmp_major = new Long64_t[fN];
Long64_t *tmp_minor = new Long64_t[fN];
Long64_t i;
Long64_t oldEntry = fTree->GetReadEntry();
Int_t current = -1;
for (i=0;i<fN;i++) {
Long64_t centry = fTree->LoadTree(i);
if (centry < 0) break;
if (fTree->GetTreeNumber() != current) {
current = fTree->GetTreeNumber();
fMajorFormula->UpdateFormulaLeaves();
fMinorFormula->UpdateFormulaLeaves();
}
tmp_major[i] = (Long64_t) fMajorFormula->EvalInstance<LongDouble_t>();
tmp_minor[i] = (Long64_t) fMinorFormula->EvalInstance<LongDouble_t>();
}
fIndex = new Long64_t[fN];
for(i = 0; i < fN; i++) { fIndex[i] = i; }
std::sort(fIndex, fIndex + fN, IndexSortComparator(tmp_major, tmp_minor) );
fIndexValues = new Long64_t[fN];
fIndexValuesMinor = new Long64_t[fN];
for (i=0;i<fN;i++) {
fIndexValues[i] = tmp_major[fIndex[i]];
fIndexValuesMinor[i] = tmp_minor[fIndex[i]];
}
delete [] tmp_major;
delete [] tmp_minor;
fTree->LoadTree(oldEntry);
}
TTreeIndex::~TTreeIndex()
{
if (fTree && fTree->GetTreeIndex() == this) fTree->SetTreeIndex(0);
delete [] fIndexValues; fIndexValues = 0;
delete [] fIndexValuesMinor; fIndexValuesMinor = 0;
delete [] fIndex; fIndex = 0;
delete fMajorFormula; fMajorFormula = 0;
delete fMinorFormula; fMinorFormula = 0;
delete fMajorFormulaParent; fMajorFormulaParent = 0;
delete fMinorFormulaParent; fMinorFormulaParent = 0;
}
void TTreeIndex::Append(const TVirtualIndex *add, Bool_t delaySort )
{
if (add && add->GetN()) {
const TTreeIndex *ti_add = dynamic_cast<const TTreeIndex*>(add);
if (ti_add == 0) {
Error("Append","Can only Append a TTreeIndex to a TTreeIndex but got a %s",
add->IsA()->GetName());
}
Long64_t oldn = fN;
fN += add->GetN();
Long64_t *oldIndex = fIndex;
Long64_t *oldValues = GetIndexValues();
Long64_t *oldValues2 = GetIndexValuesMinor();
fIndex = new Long64_t[fN];
fIndexValues = new Long64_t[fN];
fIndexValuesMinor = new Long64_t[fN];
Long_t size = sizeof(Long64_t) * oldn;
Long_t add_size = sizeof(Long64_t) * add->GetN();
memcpy(fIndex,oldIndex, size);
memcpy(fIndexValues,oldValues, size);
memcpy(fIndexValuesMinor,oldValues2, size);
Long64_t *addIndex = ti_add->GetIndex();
Long64_t *addValues = ti_add->GetIndexValues();
Long64_t *addValues2 = ti_add->GetIndexValuesMinor();
memcpy(fIndex + oldn, addIndex, add_size);
memcpy(fIndexValues + oldn, addValues, add_size);
memcpy(fIndexValuesMinor + oldn, addValues2, add_size);
for(Int_t i = 0; i < add->GetN(); i++) {
fIndex[oldn + i] += oldn;
}
delete [] oldIndex;
delete [] oldValues;
delete [] oldValues2;
}
if (!delaySort) {
Long64_t *addValues = GetIndexValues();
Long64_t *addValues2 = GetIndexValuesMinor();
Long64_t *ind = fIndex;
Long64_t *conv = new Long64_t[fN];
for(Long64_t i = 0; i < fN; i++) { conv[i] = i; }
std::sort(conv, conv+fN, IndexSortComparator(addValues, addValues2) );
fIndex = new Long64_t[fN];
fIndexValues = new Long64_t[fN];
fIndexValuesMinor = new Long64_t[fN];
for (Int_t i=0;i<fN;i++) {
fIndex[i] = ind[conv[i]];
fIndexValues[i] = addValues[conv[i]];
fIndexValuesMinor[i] = addValues2[conv[i]];
}
delete [] addValues;
delete [] addValues2;
delete [] ind;
delete [] conv;
}
}
bool TTreeIndex::ConvertOldToNew()
{
if( !fIndexValuesMinor && fN ) {
fIndexValuesMinor = new Long64_t[fN];
for(int i=0; i<fN; i++) {
fIndexValuesMinor[i] = (fIndexValues[i] & 0x7fffffff);
fIndexValues[i] >>= 31;
}
return true;
}
return false;
}
Long64_t TTreeIndex::GetEntryNumberFriend(const TTree *parent)
{
if (!parent) return -3;
GetMajorFormulaParent(parent);
GetMinorFormulaParent(parent);
if (!fMajorFormulaParent || !fMinorFormulaParent) return -1;
if (!fMajorFormulaParent->GetNdim() || !fMinorFormulaParent->GetNdim()) {
Long64_t pentry = parent->GetReadEntry();
if (pentry >= fTree->GetEntries()) return -2;
return pentry;
}
Double_t majord = fMajorFormulaParent->EvalInstance();
Double_t minord = fMinorFormulaParent->EvalInstance();
Long64_t majorv = (Long64_t)majord;
Long64_t minorv = (Long64_t)minord;
return fTree->GetEntryNumberWithIndex(majorv,minorv);
}
Long64_t TTreeIndex::FindValues(Long64_t major, Long64_t minor) const
{
Long64_t mid, step, pos = 0, count = fN;
while( count > 0 ) {
step = count / 2;
mid = pos + step;
if( fIndexValues[mid] < major
|| ( fIndexValues[mid] == major && fIndexValuesMinor[mid] < minor ) ) {
pos = mid+1;
count -= step + 1;
} else
count = step;
}
return pos;
}
Long64_t TTreeIndex::GetEntryNumberWithBestIndex(Long64_t major, Long64_t minor) const
{
if (fN == 0) return -1;
Long64_t pos = FindValues(major, minor);
if( pos < fN && fIndexValues[pos] == major && fIndexValuesMinor[pos] == minor )
return fIndex[pos];
if( --pos < 0 )
return -1;
return fIndex[pos];
}
Long64_t TTreeIndex::GetEntryNumberWithIndex(Long64_t major, Long64_t minor) const
{
if (fN == 0) return -1;
Long64_t pos = FindValues(major, minor);
if( pos < fN && fIndexValues[pos] == major && fIndexValuesMinor[pos] == minor )
return fIndex[pos];
return -1;
}
Long64_t* TTreeIndex::GetIndexValuesMinor() const
{
return fIndexValuesMinor;
}
TTreeFormula *TTreeIndex::GetMajorFormula()
{
if (!fMajorFormula) {
fMajorFormula = new TTreeFormula("Major",fMajorName.Data(),fTree);
fMajorFormula->SetQuickLoad(kTRUE);
}
return fMajorFormula;
}
TTreeFormula *TTreeIndex::GetMinorFormula()
{
if (!fMinorFormula) {
fMinorFormula = new TTreeFormula("Minor",fMinorName.Data(),fTree);
fMinorFormula->SetQuickLoad(kTRUE);
}
return fMinorFormula;
}
TTreeFormula *TTreeIndex::GetMajorFormulaParent(const TTree *parent)
{
if (!fMajorFormulaParent) {
TTree::TFriendLock friendlock(fTree, TTree::kFindLeaf | TTree::kFindBranch | TTree::kGetBranch | TTree::kGetLeaf);
fMajorFormulaParent = new TTreeFormula("MajorP",fMajorName.Data(),const_cast<TTree*>(parent));
fMajorFormulaParent->SetQuickLoad(kTRUE);
}
if (fMajorFormulaParent->GetTree() != parent) {
fMajorFormulaParent->SetTree(const_cast<TTree*>(parent));
fMajorFormulaParent->UpdateFormulaLeaves();
}
return fMajorFormulaParent;
}
TTreeFormula *TTreeIndex::GetMinorFormulaParent(const TTree *parent)
{
if (!fMinorFormulaParent) {
TTree::TFriendLock friendlock(fTree, TTree::kFindLeaf | TTree::kFindBranch | TTree::kGetBranch | TTree::kGetLeaf);
fMinorFormulaParent = new TTreeFormula("MinorP",fMinorName.Data(),const_cast<TTree*>(parent));
fMinorFormulaParent->SetQuickLoad(kTRUE);
}
if (fMinorFormulaParent->GetTree() != parent) {
fMinorFormulaParent->SetTree(const_cast<TTree*>(parent));
fMinorFormulaParent->UpdateFormulaLeaves();
}
return fMinorFormulaParent;
}
void TTreeIndex::Print(Option_t * option) const
{
TString opt = option;
Bool_t printEntry = kFALSE;
Long64_t n = fN;
if (opt.Contains("10")) n = 10;
if (opt.Contains("100")) n = 100;
if (opt.Contains("1000")) n = 1000;
if (opt.Contains("all")) {
printEntry = kTRUE;
}
if (printEntry) {
Printf("\n*****************************************************************");
Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
Printf("*****************************************************************");
Printf("%8s : %16s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data(),"entry number");
Printf("*****************************************************************");
for (Long64_t i=0;i<n;i++) {
Printf("%8lld : %8lld : %8lld : %8lld",
i, fIndexValues[i], GetIndexValuesMinor()[i], fIndex[i]);
}
} else {
Printf("\n**********************************************");
Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
Printf("**********************************************");
Printf("%8s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data());
Printf("**********************************************");
for (Long64_t i=0;i<n;i++) {
Printf("%8lld : %8lld : %8lld",
i, fIndexValues[i],GetIndexValuesMinor()[i]);
}
}
}
void TTreeIndex::Streamer(TBuffer &R__b)
{
UInt_t R__s, R__c;
if (R__b.IsReading()) {
Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
TVirtualIndex::Streamer(R__b);
fMajorName.Streamer(R__b);
fMinorName.Streamer(R__b);
R__b >> fN;
fIndexValues = new Long64_t[fN];
R__b.ReadFastArray(fIndexValues,fN);
if( R__v > 1 ) {
fIndexValuesMinor = new Long64_t[fN];
R__b.ReadFastArray(fIndexValuesMinor,fN);
} else {
ConvertOldToNew();
}
fIndex = new Long64_t[fN];
R__b.ReadFastArray(fIndex,fN);
R__b.CheckByteCount(R__s, R__c, TTreeIndex::IsA());
} else {
R__c = R__b.WriteVersion(TTreeIndex::IsA(), kTRUE);
TVirtualIndex::Streamer(R__b);
fMajorName.Streamer(R__b);
fMinorName.Streamer(R__b);
R__b << fN;
R__b.WriteFastArray(fIndexValues, fN);
R__b.WriteFastArray(fIndexValuesMinor, fN);
R__b.WriteFastArray(fIndex, fN);
R__b.SetByteCount(R__c, kTRUE);
}
}
void TTreeIndex::UpdateFormulaLeaves(const TTree *parent)
{
if (fMajorFormula) { fMajorFormula->UpdateFormulaLeaves();}
if (fMinorFormula) { fMinorFormula->UpdateFormulaLeaves();}
if (fMajorFormulaParent) {
if (parent) fMajorFormulaParent->SetTree(const_cast<TTree*>(parent));
fMajorFormulaParent->UpdateFormulaLeaves();
}
if (fMinorFormulaParent) {
if (parent) fMinorFormulaParent->SetTree(const_cast<TTree*>(parent));
fMinorFormulaParent->UpdateFormulaLeaves();
}
}
void TTreeIndex::SetTree(const TTree *T)
{
fTree = (TTree*)T;
}