#include "TChainIndex.h"
#include "TChain.h"
#include "TTreeFormula.h"
#include "TTreeIndex.h"
#include "TFile.h"
#include "TError.h"
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");
}
R__ASSERT(dynamic_cast<TTreeIndex*>(index));
entry.fMinIndexValue = dynamic_cast<TTreeIndex*>(index)->GetIndexValues()[0];
entry.fMaxIndexValue = dynamic_cast<TTreeIndex*>(index)->GetIndexValues()[index->GetN() - 1];
fEntries.push_back(entry);
}
for (i = 0; i < Int_t(fEntries.size() - 1); i++) {
if (fEntries[i].fMaxIndexValue > fEntries[i+1].fMinIndexValue) {
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) {
R__ASSERT(dynamic_cast<const TTreeIndex*>(index));
TChainIndexEntry entry;
entry.fTreeIndex = 0;
entry.fMinIndexValue = dynamic_cast<const TTreeIndex*>(index)->GetIndexValues()[0];
entry.fMaxIndexValue = dynamic_cast<const TTreeIndex*>(index)->GetIndexValues()[index->GetN() - 1];
fEntries.push_back(entry);
}
if (!delaySort) {
for (Int_t i = 0; i < Int_t(fEntries.size() - 1); i++) {
if (fEntries[i].fMaxIndexValue > fEntries[i+1].fMinIndexValue) {
DeleteIndices();
MakeZombie();
Error("TChainIndex", "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(Int_t major, Int_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);
}
Long64_t indexValue = (Long64_t(major) << 31) + minor;
if (indexValue < fEntries[0].fMinIndexValue) {
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].fMinIndexValue) {
treeNo = i;
break;
}
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);
}
}
Int_t TChainIndex::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 TChainIndex::GetEntryNumberWithBestIndex(Int_t major, Int_t minor) const
{
std::pair<TVirtualIndex*, Int_t> indexAndNumber = GetSubTreeIndex(major, minor);
if (!indexAndNumber.first) {
Error("GetEntryNumberWithBestIndex","no index found");
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(Int_t major, Int_t minor) const
{
std::pair<TVirtualIndex*, Int_t> indexAndNumber = GetSubTreeIndex(major, minor);
if (!indexAndNumber.first) {
Error("GetEntryNumberWithIndex","no index found");
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);
return rv + chain->GetTreeOffset()[indexAndNumber.second];
}
}
TTreeFormula *TChainIndex::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 *TChainIndex::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 TChainIndex::UpdateFormulaLeaves(const TTree *parent)
{
if (fMajorFormulaParent) {
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);
}
Last change: Wed Dec 10 16:17:44 2008
Last generated: 2008-12-10 16:17
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.