#include "TTreeIndex.h"
#include "TTree.h"
#include "TMath.h"
ClassImp(TTreeIndex)
TTreeIndex::TTreeIndex(): TVirtualIndex()
{
fTree = 0;
fN = 0;
fIndexValues = 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;
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 *w = 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();
}
Double_t majord = fMajorFormula->EvalInstance();
Double_t minord = fMinorFormula->EvalInstance();
Long64_t majorv = (Long64_t)majord;
Long64_t minorv = (Long64_t)minord;
w[i] = majorv<<31;
w[i] += minorv;
}
fIndex = new Long64_t[fN];
TMath::Sort(fN,w,fIndex,0);
fIndexValues = new Long64_t[fN];
for (i=0;i<fN;i++) {
fIndexValues[i] = w[fIndex[i]];
}
delete [] w;
fTree->LoadTree(oldEntry);
}
TTreeIndex::~TTreeIndex()
{
if (fTree && fTree->GetTreeIndex() == this) fTree->SetTreeIndex(0);
delete [] fIndexValues; fIndexValues = 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()) {
Long64_t oldn = fN;
fN += add->GetN();
Long64_t *oldIndex = fIndex;
Long64_t *oldValues = fIndexValues;
fIndex = new Long64_t[fN];
fIndexValues = 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);
Long64_t *addIndex = dynamic_cast<const TTreeIndex*>(add)->GetIndex();
Long64_t *addValues = dynamic_cast<const TTreeIndex*>(add)->GetIndexValues();
memcpy(fIndex + oldn, addIndex, add_size);
memcpy(fIndexValues + oldn, addValues, add_size);
for(Int_t i = 0; i < add->GetN(); i++) {
fIndex[oldn + i] += oldn;
}
delete [] oldIndex;
delete [] oldValues;
}
if (!delaySort) {
Long64_t *w = fIndexValues;
Long64_t *ind = fIndex;
Long64_t *conv = new Long64_t[fN];
TMath::Sort(fN,w,conv,0);
fIndex = new Long64_t[fN];
fIndexValues = new Long64_t[fN];
for (Int_t i=0;i<fN;i++) {
fIndex[i] = ind[conv[i]];
fIndexValues[i] = w[conv[i]];
}
delete [] w;
delete [] ind;
delete [] conv;
}
}
Int_t TTreeIndex::GetEntryNumberFriend(const TTree *T)
{
if (!T) return -3;
GetMajorFormulaParent(T);
GetMinorFormulaParent(T);
if (!fMajorFormulaParent || !fMinorFormulaParent) return -1;
if (!fMajorFormulaParent->GetNdim() || !fMinorFormulaParent->GetNdim()) {
Int_t pentry = T->GetReadEntry();
if (pentry >= (Int_t)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::GetEntryNumberWithBestIndex(Int_t major, Int_t minor) const
{
if (fN == 0) return -1;
Long64_t value = Long64_t(major)<<31;
value += minor;
Int_t i = TMath::BinarySearch(fN, fIndexValues, value);
if (i < 0) return -1;
return fIndex[i];
}
Long64_t TTreeIndex::GetEntryNumberWithIndex(Int_t major, Int_t minor) const
{
if (fN == 0) return -1;
Long64_t value = Long64_t(major)<<31;
value += minor;
Int_t i = TMath::BinarySearch(fN, fIndexValues, value);
if (i < 0) return -1;
if (fIndexValues[i] != value) return -1;
return fIndex[i];
}
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 *T)
{
if (!fMajorFormulaParent) {
fMajorFormulaParent = new TTreeFormula("MajorP",fMajorName.Data(),(TTree*)T);
fMajorFormulaParent->SetQuickLoad(kTRUE);
}
if (fMajorFormulaParent->GetTree() != T) {
fMajorFormulaParent->SetTree((TTree*)T);
fMajorFormulaParent->UpdateFormulaLeaves();
}
return fMajorFormulaParent;
}
TTreeFormula *TTreeIndex::GetMinorFormulaParent(const TTree *T)
{
if (!fMinorFormulaParent) {
fMinorFormulaParent = new TTreeFormula("MinorP",fMinorName.Data(),(TTree*)T);
fMinorFormulaParent->SetQuickLoad(kTRUE);
}
if (fMinorFormulaParent->GetTree() != T) {
fMinorFormulaParent->SetTree((TTree*)T);
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++) {
Long64_t minor = fIndexValues[i] & 0xffff;
Long64_t major = fIndexValues[i]>>31;
Printf("%8lld : %8lld : %8lld : %8lld",i,major,minor,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++) {
Long64_t minor = fIndexValues[i] & 0xffff;
Long64_t major = fIndexValues[i]>>31;
Printf("%8lld : %8lld : %8lld",i,major,minor);
}
}
}
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);
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(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((TTree*)parent);
fMajorFormulaParent->UpdateFormulaLeaves();
}
if (fMinorFormulaParent) {
if (parent) fMinorFormulaParent->SetTree((TTree*)parent);
fMinorFormulaParent->UpdateFormulaLeaves();
}
}
void TTreeIndex::SetTree(const TTree *T)
{
fTree = (TTree*)T;
}