#include "TChainIndex.h"
#include "TChain.h"
#include "TTreeFormula.h"
#include "TTreeIndex.h"
#include "TFile.h"
#include "TError.h"
void TChainIndex::TChainIndexEntry::SetMinMaxFrom(const TTreeIndex *index )
{
fMinIndexValue = index->GetIndexValues()[0];
fMinIndexValMinor = index->GetIndexValuesMinor()[0];
fMaxIndexValue = index->GetIndexValues()[index->GetN() - 1];
fMaxIndexValMinor = index->GetIndexValuesMinor()[index->GetN() - 1];
}
ClassImp(TChainIndex)
TChainIndex::TChainIndex(): TVirtualIndex()
{
fTree = 0;
fMajorFormulaParent = fMinorFormulaParent = 0;
}
TChainIndex::TChainIndex(const TTree *T, const char *majorname, const char *minorname)
: TVirtualIndex()
{
fTree = 0;
fMajorFormulaParent = fMinorFormulaParent = 0;
TChain *chain = dynamic_cast<TChain*>(const_cast<TTree*>(T));
if (!chain) {
MakeZombie();
Error("TChainIndex", "Cannot create a TChainIndex."
" The Tree passed as an argument is not a TChain");
return;
}
fTree = (TTree*)T;
fMajorName = majorname;
fMinorName = minorname;
Int_t i = 0;
for (i = 0; i < chain->GetNtrees(); i++) {
chain->LoadTree((chain->GetTreeOffset())[i]);
TVirtualIndex *index = chain->GetTree()->GetTreeIndex();
TChainIndexEntry entry;
entry.fTreeIndex = 0;
if (index) {
if (strcmp(majorname,index->GetMajorName()) || strcmp(minorname,index->GetMinorName())) {
MakeZombie();
Error("TChainIndex","Tree in file %s has an index built with majorname=%s and minorname=%s",chain->GetTree()->GetCurrentFile()->GetName(),index->GetMajorName(),index->GetMinorName());
return;
}
}
if (!index) {
chain->GetTree()->BuildIndex(majorname, minorname);
index = chain->GetTree()->GetTreeIndex();
chain->GetTree()->SetTreeIndex(0);
entry.fTreeIndex = index;
}
if (!index || index->IsZombie() || index->GetN() == 0) {
DeleteIndices();
MakeZombie();
Error("TChainIndex", "Error creating a tree index on a tree in the chain");
return;
}
TTreeIndex *ti_index = dynamic_cast<TTreeIndex*>(index);
if (ti_index == 0) {
Error("TChainIndex", "The underlying TTree must have a TTreeIndex but has a %s.",
index->IsA()->GetName());
return;
}
entry.SetMinMaxFrom(ti_index);
fEntries.push_back(entry);
}
for (i = 0; i < Int_t(fEntries.size() - 1); i++) {
if( fEntries[i].GetMaxIndexValPair() > fEntries[i+1].GetMinIndexValPair() ) {
DeleteIndices();
MakeZombie();
Error("TChainIndex", "The indices in files of this chain aren't sorted.");
}
}
}
void TChainIndex::Append(const TVirtualIndex *index, Bool_t delaySort )
{
if (index) {
const TTreeIndex *ti_index = dynamic_cast<const TTreeIndex*>(index);
if (ti_index == 0) {
Error("Append", "The given index is not a TTreeIndex but a %s",
index->IsA()->GetName());
}
TChainIndexEntry entry;
entry.fTreeIndex = 0;
entry.SetMinMaxFrom(ti_index);
fEntries.push_back(entry);
}
if (!delaySort) {
for (Int_t i = 0; i < Int_t(fEntries.size() - 1); i++) {
if( fEntries[i].GetMaxIndexValPair() > fEntries[i+1].GetMinIndexValPair() ) {
DeleteIndices();
MakeZombie();
Error("Append", "The indices in files of this chain aren't sorted.");
}
}
}
}
void TChainIndex::DeleteIndices()
{
for (unsigned int i = 0; i < fEntries.size(); i++) {
if (fEntries[i].fTreeIndex) {
if (fTree->GetTree() && fTree->GetTree()->GetTreeIndex() == fEntries[i].fTreeIndex) {
fTree->GetTree()->SetTreeIndex(0);
SafeDelete(fEntries[i].fTreeIndex);
}
SafeDelete(fEntries[i].fTreeIndex);
}
}
}
TChainIndex::~TChainIndex()
{
DeleteIndices();
if (fTree && fTree->GetTreeIndex() == this)
fTree->SetTreeIndex(0);
}
std::pair<TVirtualIndex*, Int_t> TChainIndex::GetSubTreeIndex(Long64_t major, Long64_t minor) const
{
using namespace std;
if (fEntries.size() == 0) {
Warning("GetSubTreeIndex", "No subindices in the chain. The chain is probably empty");
return make_pair(static_cast<TVirtualIndex*>(0), 0);
}
const TChainIndexEntry::IndexValPair_t indexValue(major, minor);
if( indexValue < fEntries[0].GetMinIndexValPair() ) {
Warning("GetSubTreeIndex", "The index value is less than the smallest index values in subtrees");
return make_pair(static_cast<TVirtualIndex*>(0), 0);
}
Int_t treeNo = fEntries.size() - 1;
for (unsigned int i = 0; i < fEntries.size() - 1; i++) {
if( indexValue < fEntries[i+1].GetMinIndexValPair() ) {
treeNo = i;
break;
}
}
if( indexValue > fEntries[treeNo].GetMaxIndexValPair() ) {
return make_pair(static_cast<TVirtualIndex*>(0), 0);
}
TChain* chain = dynamic_cast<TChain*> (fTree);
R__ASSERT(chain);
chain->LoadTree(chain->GetTreeOffset()[treeNo]);
TVirtualIndex* index = fTree->GetTree()->GetTreeIndex();
if (index)
return make_pair(static_cast<TVirtualIndex*>(index), treeNo);
else {
index = fEntries[treeNo].fTreeIndex;
if (!index) {
Warning("GetSubTreeIndex", "The tree has no index and the chain index"
" doesn't store an index for that tree");
return make_pair(static_cast<TVirtualIndex*>(0), 0);
}
else {
fTree->GetTree()->SetTreeIndex(index);
return make_pair(static_cast<TVirtualIndex*>(index), treeNo);
}
}
}
void TChainIndex::ReleaseSubTreeIndex(TVirtualIndex* index, int treeNo) const
{
if (fEntries[treeNo].fTreeIndex == index) {
R__ASSERT(fTree->GetTree()->GetTreeIndex() == index);
fTree->GetTree()->SetTreeIndex(0);
}
}
Long64_t TChainIndex::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 TChainIndex::GetEntryNumberWithBestIndex(Long64_t major, Long64_t minor) const
{
std::pair<TVirtualIndex*, Int_t> indexAndNumber = GetSubTreeIndex(major, minor);
if (!indexAndNumber.first) {
return -1;
}
else {
Long64_t rv = indexAndNumber.first->GetEntryNumberWithBestIndex(major, minor);
ReleaseSubTreeIndex(indexAndNumber.first, indexAndNumber.second);
TChain* chain = dynamic_cast<TChain*> (fTree);
R__ASSERT(chain);
return rv + chain->GetTreeOffset()[indexAndNumber.second];
}
}
Long64_t TChainIndex::GetEntryNumberWithIndex(Long64_t major, Long64_t minor) const
{
std::pair<TVirtualIndex*, Int_t> indexAndNumber = GetSubTreeIndex(major, minor);
if (!indexAndNumber.first) {
return -1;
}
else {
Long64_t rv = indexAndNumber.first->GetEntryNumberWithIndex(major, minor);
ReleaseSubTreeIndex(indexAndNumber.first, indexAndNumber.second);
TChain* chain = dynamic_cast<TChain*> (fTree);
R__ASSERT(chain);
if (rv >= 0) {
return rv + chain->GetTreeOffset()[indexAndNumber.second];
} else {
return rv;
}
}
}
TTreeFormula *TChainIndex::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 *TChainIndex::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 TChainIndex::UpdateFormulaLeaves(const TTree *parent)
{
if (fMajorFormulaParent) {
TTree::TFriendLock friendlock(fTree, TTree::kFindLeaf | TTree::kFindBranch | TTree::kGetBranch | TTree::kGetLeaf);
if (parent) fMajorFormulaParent->SetTree((TTree*)parent);
fMajorFormulaParent->UpdateFormulaLeaves();
}
if (fMinorFormulaParent) {
if (parent) fMinorFormulaParent->SetTree((TTree*)parent);
fMinorFormulaParent->UpdateFormulaLeaves();
}
}
void TChainIndex::SetTree(const TTree *T)
{
R__ASSERT(fTree == 0 || fTree == T || T==0);
}