TTreeIndex.cxx

Go to the documentation of this file.
00001 // @(#)root/tree:$Id: TTreeIndex.cxx 35344 2010-09-16 21:34:21Z pcanal $
00002 // Author: Rene Brun   05/07/2004
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //////////////////////////////////////////////////////////////////////////
00013 //                                                                      //
00014 // A Tree Index with majorname and minorname.                           //
00015 //                                                                      //
00016 //////////////////////////////////////////////////////////////////////////
00017 
00018 #include "TTreeIndex.h"
00019 #include "TTree.h"
00020 #include "TMath.h"
00021 
00022 ClassImp(TTreeIndex)
00023 
00024 //______________________________________________________________________________
00025 TTreeIndex::TTreeIndex(): TVirtualIndex()
00026 {
00027    // Default constructor for TTreeIndex
00028 
00029    fTree               = 0;
00030    fN                  = 0;
00031    fIndexValues        = 0;
00032    fIndex              = 0;
00033    fMajorFormula       = 0;
00034    fMinorFormula       = 0;
00035    fMajorFormulaParent = 0;
00036    fMinorFormulaParent = 0;
00037 }
00038 
00039 //______________________________________________________________________________
00040 TTreeIndex::TTreeIndex(const TTree *T, const char *majorname, const char *minorname)
00041            : TVirtualIndex()
00042 {
00043    // Normal constructor for TTreeIndex
00044    //
00045    // Build an index table using the leaves of Tree T with  major & minor names
00046    // The index is built with the expressions given in "majorname" and "minorname".
00047    //
00048    // a Long64_t array fIndexValues is built with:
00049    //    major = the value of majorname converted to an integer
00050    //    minor = the value of minorname converted to an integer
00051    //    fIndexValues[i] = major<<31 + minor
00052    // This array is sorted. The sorted fIndex[i] contains the serial number
00053    // in the Tree corresponding to the pair "major,minor" in fIndexvalues[i].
00054    //
00055    //  Once the index is computed, one can retrieve one entry via
00056    //    T->GetEntryWithIndex(majornumber, minornumber)
00057    // Example:
00058    //  tree.BuildIndex("Run","Event"); //creates an index using leaves Run and Event
00059    //  tree.GetEntryWithIndex(1234,56789); //reads entry corresponding to
00060    //                                        Run=1234 and Event=56789
00061    //
00062    // Note that majorname and minorname may be expressions using original
00063    // Tree variables eg: "run-90000", "event +3*xx". However the result
00064    // must be integer.
00065    // In case an expression is specified, the equivalent expression must be computed
00066    // when calling GetEntryWithIndex.
00067    //
00068    // To build an index with only majorname, specify minorname="0" (default)
00069    //
00070    //    TreeIndex and Friend Trees
00071    //    ---------------------------
00072    // Assuming a parent Tree T and a friend Tree TF, the following cases are supported:
00073    // CASE 1: T->GetEntry(entry) is called
00074    //         In this case, the serial number entry is used to retrieve
00075    //         the data in both Trees.
00076    // CASE 2: T->GetEntry(entry) is called, TF has a TreeIndex
00077    //         the expressions given in major/minorname of TF are used
00078    //         to compute the value pair major,minor with the data in T.
00079    //         TF->GetEntryWithIndex(major,minor) is then called (tricky case!)
00080    // CASE 3: T->GetEntryWithIndex(major,minor) is called.
00081    //         It is assumed that both T and TF have a TreeIndex built using
00082    //         the same major and minor name.
00083    //
00084    //    Saving the TreeIndex
00085    //    --------------------
00086    // Once the index is built, it can be saved with the TTree object
00087    // with tree.Write(); (if the file has been open in "update" mode).
00088    //
00089    // The most convenient place to create the index is at the end of
00090    // the filling process just before saving the Tree header.
00091    // If a previous index was computed, it is redefined by this new call.
00092    //
00093    // Note that this function can also be applied to a TChain.
00094    //
00095    // The return value is the number of entries in the Index (< 0 indicates failure)
00096    //
00097    // It is possible to play with different TreeIndex in the same Tree.
00098    // see comments in TTree::SetTreeIndex.
00099 
00100    fTree               = (TTree*)T;
00101    fN                  = 0;
00102    fIndexValues        = 0;
00103    fIndex              = 0;
00104    fMajorFormula       = 0;
00105    fMinorFormula       = 0;
00106    fMajorFormulaParent = 0;
00107    fMinorFormulaParent = 0;
00108    fMajorName          = majorname;
00109    fMinorName          = minorname;
00110    if (!T) return;
00111    fN = T->GetEntries();
00112    if (fN <= 0) {
00113       MakeZombie();
00114       Error("TreeIndex","Cannot build a TreeIndex with a Tree having no entries");
00115       return;
00116    }
00117 
00118    GetMajorFormula();
00119    GetMinorFormula();
00120    if (!fMajorFormula || !fMinorFormula) {
00121       MakeZombie();
00122       Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
00123       return;
00124    }
00125    if ((fMajorFormula->GetNdim() != 1) || (fMinorFormula->GetNdim() != 1)) {
00126       MakeZombie();
00127       Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
00128       return;
00129    }
00130    // accessing array elements should be OK
00131    //if ((fMajorFormula->GetMultiplicity() != 0) || (fMinorFormula->GetMultiplicity() != 0)) {
00132    //   MakeZombie();
00133    //   Error("TreeIndex","Cannot build the index with major=%s, minor=%s that cannot be arrays",fMajorName.Data(), fMinorName.Data());
00134    //   return;
00135    //}
00136 
00137    Long64_t *w = new Long64_t[fN];
00138    Long64_t i;
00139    Long64_t oldEntry = fTree->GetReadEntry();
00140    Int_t current = -1;
00141    for (i=0;i<fN;i++) {
00142       Long64_t centry = fTree->LoadTree(i);
00143       if (centry < 0) break;
00144       if (fTree->GetTreeNumber() != current) {
00145          current = fTree->GetTreeNumber();
00146          fMajorFormula->UpdateFormulaLeaves();
00147          fMinorFormula->UpdateFormulaLeaves();
00148       }
00149       Double_t majord = fMajorFormula->EvalInstance();
00150       Double_t minord = fMinorFormula->EvalInstance();
00151       Long64_t majorv = (Long64_t)majord;
00152       Long64_t minorv = (Long64_t)minord;
00153       w[i]  = majorv<<31;
00154       w[i] += minorv;
00155    }
00156    fIndex = new Long64_t[fN];
00157    TMath::Sort(fN,w,fIndex,0);
00158    fIndexValues = new Long64_t[fN];
00159    for (i=0;i<fN;i++) {
00160       fIndexValues[i] = w[fIndex[i]];
00161    }
00162 
00163    delete [] w;
00164    fTree->LoadTree(oldEntry);
00165 }
00166 
00167 //______________________________________________________________________________
00168 TTreeIndex::~TTreeIndex()
00169 {
00170    // Destructor.
00171 
00172    if (fTree && fTree->GetTreeIndex() == this) fTree->SetTreeIndex(0);
00173    delete [] fIndexValues;      fIndexValues = 0;
00174    delete [] fIndex;            fIndex = 0;
00175    delete fMajorFormula;        fMajorFormula  = 0;
00176    delete fMinorFormula;        fMinorFormula  = 0;
00177    delete fMajorFormulaParent;  fMajorFormulaParent = 0;
00178    delete fMinorFormulaParent;  fMinorFormulaParent = 0;
00179 }
00180 
00181 //______________________________________________________________________________
00182 void TTreeIndex::Append(const TVirtualIndex *add, Bool_t delaySort ) 
00183 {
00184    // Append 'add' to this index.  Entry 0 in add will become entry n+1 in this.
00185    // If delaySort is true, do not sort the value, then you must call
00186    // Append(0,kFALSE);
00187 
00188    
00189    if (add && add->GetN()) {
00190       // Create new buffer (if needed)
00191 
00192       const TTreeIndex *ti_add = dynamic_cast<const TTreeIndex*>(add);
00193       if (ti_add == 0) {
00194          Error("Append","Can only Append a TTreeIndex to a TTreeIndex but got a %s",
00195                add->IsA()->GetName());
00196       }
00197 
00198       Long64_t oldn = fN;
00199       fN += add->GetN();
00200 
00201       Long64_t *oldIndex = fIndex;
00202       Long64_t *oldValues = fIndexValues;
00203 
00204       fIndex = new Long64_t[fN];
00205       fIndexValues = new Long64_t[fN];
00206 
00207       // Copy data
00208 
00209       Long_t size = sizeof(Long64_t) * oldn;
00210       Long_t add_size = sizeof(Long64_t) * add->GetN();
00211       
00212       memcpy(fIndex,oldIndex, size);
00213       memcpy(fIndexValues,oldValues, size);
00214 
00215       Long64_t *addIndex = ti_add->GetIndex();
00216       Long64_t *addValues = ti_add->GetIndexValues();
00217 
00218       memcpy(fIndex + oldn, addIndex, add_size);
00219       memcpy(fIndexValues + oldn, addValues, add_size);
00220       for(Int_t i = 0; i < add->GetN(); i++) {
00221          fIndex[oldn + i] += oldn;
00222       }
00223 
00224       delete [] oldIndex;
00225       delete [] oldValues;
00226    }
00227 
00228    // Sort.
00229    if (!delaySort) {
00230       Long64_t *w = fIndexValues;
00231       Long64_t *ind = fIndex;
00232       Long64_t *conv = new Long64_t[fN];
00233 
00234       TMath::Sort(fN,w,conv,0);
00235 
00236       fIndex = new Long64_t[fN];
00237       fIndexValues = new Long64_t[fN];
00238 
00239       for (Int_t i=0;i<fN;i++) {
00240          fIndex[i] = ind[conv[i]];
00241          fIndexValues[i] = w[conv[i]];
00242       }
00243       delete [] w;
00244       delete [] ind;
00245       delete [] conv;
00246    }
00247 }
00248 
00249 //______________________________________________________________________________
00250 Int_t TTreeIndex::GetEntryNumberFriend(const TTree *T)
00251 {
00252 // returns the entry number in this friend Tree corresponding to entry in
00253 // the master Tree T.
00254 // In case this friend Tree and T do not share an index with the same
00255 // major and minor name, the entry serial number in the friend tree
00256 // and in the master Tree are assumed to be the same
00257 
00258    if (!T) return -3;
00259    GetMajorFormulaParent(T);
00260    GetMinorFormulaParent(T);
00261    if (!fMajorFormulaParent || !fMinorFormulaParent) return -1;
00262    if (!fMajorFormulaParent->GetNdim() || !fMinorFormulaParent->GetNdim()) {
00263       // The Tree Index in the friend has a pair majorname,minorname
00264       // not available in the parent Tree T.
00265       // if the friend Tree has less entries than the parent, this is an error
00266       Int_t pentry = T->GetReadEntry();
00267       if (pentry >= (Int_t)fTree->GetEntries()) return -2;
00268       // otherwise we ignore the Tree Index and return the entry number
00269       // in the parent Tree.
00270       return pentry;
00271    }
00272 
00273    // majorname, minorname exist in the parent Tree
00274    // we find the current values pair majorv,minorv in the parent Tree
00275    Double_t majord = fMajorFormulaParent->EvalInstance();
00276    Double_t minord = fMinorFormulaParent->EvalInstance();
00277    Long64_t majorv = (Long64_t)majord;
00278    Long64_t minorv = (Long64_t)minord;
00279    // we check if this pair exist in the index.
00280    // if yes, we return the corresponding entry number
00281    // if not the function returns -1
00282    return fTree->GetEntryNumberWithIndex(majorv,minorv);
00283 }
00284 
00285 //______________________________________________________________________________
00286 Long64_t TTreeIndex::GetEntryNumberWithBestIndex(Int_t major, Int_t minor) const
00287 {
00288 // Return entry number corresponding to major and minor number
00289 // Note that this function returns only the entry number, not the data
00290 // To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
00291 // the BuildIndex function has created a table of Double_t* of sorted values
00292 // corresponding to val = major<<31 + minor;
00293 // The function performs binary search in this sorted table.
00294 // If it finds a pair that maches val, it returns directly the
00295 // index in the table.
00296 // If an entry corresponding to major and minor is not found, the function
00297 // returns the index of the major,minor pair immediatly lower than the
00298 // requested value, ie it will return -1 if the pair is lower than
00299 // the first entry in the index.
00300 //
00301 // See also GetEntryNumberWithIndex
00302 
00303    if (fN == 0) return -1;
00304    Long64_t value = Long64_t(major)<<31;
00305    value += minor;
00306    Int_t i = TMath::BinarySearch(fN, fIndexValues, value);
00307    if (i < 0) return -1;
00308    return fIndex[i];
00309 }
00310 
00311 
00312 //______________________________________________________________________________
00313 Long64_t TTreeIndex::GetEntryNumberWithIndex(Int_t major, Int_t minor) const
00314 {
00315 // Return entry number corresponding to major and minor number
00316 // Note that this function returns only the entry number, not the data
00317 // To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
00318 // the BuildIndex function has created a table of Double_t* of sorted values
00319 // corresponding to val = major<<31 + minor;
00320 // The function performs binary search in this sorted table.
00321 // If it finds a pair that maches val, it returns directly the
00322 // index in the table, otherwise it returns -1.
00323 //
00324 // See also GetEntryNumberWithBestIndex
00325 
00326    if (fN == 0) return -1;
00327    Long64_t value = Long64_t(major)<<31;
00328    value += minor;
00329    Int_t i = TMath::BinarySearch(fN, fIndexValues, value);
00330    if (i < 0) return -1;
00331    if (fIndexValues[i] != value) return -1;
00332    return fIndex[i];
00333 }
00334 
00335 //______________________________________________________________________________
00336 TTreeFormula *TTreeIndex::GetMajorFormula()
00337 {
00338    // return a pointer to the TreeFormula corresponding to the majorname
00339 
00340    if (!fMajorFormula) {
00341       fMajorFormula = new TTreeFormula("Major",fMajorName.Data(),fTree);
00342       fMajorFormula->SetQuickLoad(kTRUE);
00343    }
00344    return fMajorFormula;
00345 }
00346 
00347 //______________________________________________________________________________
00348 TTreeFormula *TTreeIndex::GetMinorFormula()
00349 {
00350    // return a pointer to the TreeFormula corresponding to the minorname
00351 
00352    if (!fMinorFormula) {
00353       fMinorFormula = new TTreeFormula("Minor",fMinorName.Data(),fTree);
00354       fMinorFormula->SetQuickLoad(kTRUE);
00355    }
00356    return fMinorFormula;
00357 }
00358 
00359 //______________________________________________________________________________
00360 TTreeFormula *TTreeIndex::GetMajorFormulaParent(const TTree *T)
00361 {
00362    // return a pointer to the TreeFormula corresponding to the majorname in parent tree T
00363 
00364    if (!fMajorFormulaParent) {
00365       fMajorFormulaParent = new TTreeFormula("MajorP",fMajorName.Data(),(TTree*)T);
00366       fMajorFormulaParent->SetQuickLoad(kTRUE);
00367    }
00368    if (fMajorFormulaParent->GetTree() != T) {
00369       fMajorFormulaParent->SetTree((TTree*)T);
00370       fMajorFormulaParent->UpdateFormulaLeaves();
00371    }
00372    return fMajorFormulaParent;
00373 }
00374 
00375 //______________________________________________________________________________
00376 TTreeFormula *TTreeIndex::GetMinorFormulaParent(const TTree *T)
00377 {
00378    // return a pointer to the TreeFormula corresponding to the minorname in parent tree T
00379 
00380    if (!fMinorFormulaParent) {
00381       fMinorFormulaParent = new TTreeFormula("MinorP",fMinorName.Data(),(TTree*)T);
00382       fMinorFormulaParent->SetQuickLoad(kTRUE);
00383    }
00384    if (fMinorFormulaParent->GetTree() != T) {
00385       fMinorFormulaParent->SetTree((TTree*)T);
00386       fMinorFormulaParent->UpdateFormulaLeaves();
00387    }
00388 
00389    return fMinorFormulaParent;
00390 }
00391 
00392 //______________________________________________________________________________
00393 void TTreeIndex::Print(Option_t * option) const
00394 {
00395    // print the table with : serial number, majorname, minorname
00396    // if option = "10" print only the first 10 entries
00397    // if option = "100" print only the first 100 entries
00398    // if option = "1000" print only the first 1000 entries
00399 
00400    TString opt = option;
00401    Bool_t printEntry = kFALSE;
00402    Long64_t n = fN;
00403    if (opt.Contains("10"))   n = 10;
00404    if (opt.Contains("100"))  n = 100;
00405    if (opt.Contains("1000")) n = 1000;
00406    if (opt.Contains("all")) {
00407       printEntry = kTRUE;
00408    }
00409 
00410    if (printEntry) {
00411 
00412       Printf("\n*****************************************************************");
00413       Printf("*    Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
00414       Printf("*****************************************************************");
00415       Printf("%8s : %16s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data(),"entry number");
00416       Printf("*****************************************************************");
00417       for (Long64_t i=0;i<n;i++) {
00418          Long64_t minor = fIndexValues[i] & 0xffff;
00419          Long64_t major = fIndexValues[i]>>31;
00420          Printf("%8lld :         %8lld :         %8lld :         %8lld",i,major,minor,fIndex[i]);
00421       }
00422 
00423    } else {
00424 
00425       Printf("\n**********************************************");
00426       Printf("*    Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
00427       Printf("**********************************************");
00428       Printf("%8s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data());
00429       Printf("**********************************************");
00430       for (Long64_t i=0;i<n;i++) {
00431          Long64_t minor = fIndexValues[i] & 0xffff;
00432          Long64_t major = fIndexValues[i]>>31;
00433          Printf("%8lld :         %8lld :         %8lld",i,major,minor);
00434       }
00435    }
00436 }
00437 
00438 //______________________________________________________________________________
00439 void TTreeIndex::Streamer(TBuffer &R__b)
00440 {
00441    // Stream an object of class TTreeIndex.
00442    // Note that this Streamer should be changed to an automatic Streamer
00443    // once TStreamerInfo supports an index of type Long64_t
00444 
00445    UInt_t R__s, R__c;
00446    if (R__b.IsReading()) {
00447       Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
00448       TVirtualIndex::Streamer(R__b);
00449       fMajorName.Streamer(R__b);
00450       fMinorName.Streamer(R__b);
00451       R__b >> fN;
00452       fIndexValues = new Long64_t[fN];
00453       R__b.ReadFastArray(fIndexValues,fN);
00454       fIndex      = new Long64_t[fN];
00455       R__b.ReadFastArray(fIndex,fN);
00456       R__b.CheckByteCount(R__s, R__c, TTreeIndex::IsA());
00457    } else {
00458       R__c = R__b.WriteVersion(TTreeIndex::IsA(), kTRUE);
00459       TVirtualIndex::Streamer(R__b);
00460       fMajorName.Streamer(R__b);
00461       fMinorName.Streamer(R__b);
00462       R__b << fN;
00463       R__b.WriteFastArray(fIndexValues, fN);
00464       R__b.WriteFastArray(fIndex, fN);
00465       R__b.SetByteCount(R__c, kTRUE);
00466    }
00467 }
00468 
00469 //______________________________________________________________________________
00470 void TTreeIndex::UpdateFormulaLeaves(const TTree *parent)
00471 {
00472    // Called by TChain::LoadTree when the parent chain changes it's tree.
00473 
00474    if (fMajorFormula)       { fMajorFormula->UpdateFormulaLeaves();}
00475    if (fMinorFormula)       { fMinorFormula->UpdateFormulaLeaves();}
00476    if (fMajorFormulaParent) { 
00477       if (parent) fMajorFormulaParent->SetTree((TTree*)parent);
00478       fMajorFormulaParent->UpdateFormulaLeaves();
00479    }
00480    if (fMinorFormulaParent) { 
00481       if (parent) fMinorFormulaParent->SetTree((TTree*)parent);
00482       fMinorFormulaParent->UpdateFormulaLeaves();
00483    }
00484 }
00485 //______________________________________________________________________________
00486 void TTreeIndex::SetTree(const TTree *T)
00487 {
00488    // this function is called by TChain::LoadTree and TTreePlayer::UpdateFormulaLeaves
00489    // when a new Tree is loaded.
00490    // Because Trees in a TChain may have a different list of leaves, one
00491    // must update the leaves numbers in the TTreeFormula used by the TreeIndex.
00492 
00493    fTree = (TTree*)T;
00494 }
00495 

Generated on Tue Jul 5 15:43:02 2011 for ROOT_528-00b_version by  doxygen 1.5.1