Logo ROOT  
Reference Guide
TTreeIndex.cxx
Go to the documentation of this file.
1 // @(#)root/tree:$Id$
2 // Author: Rene Brun 05/07/2004
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class TTreeIndex
13 A Tree Index with majorname and minorname.
14 */
15 
16 #include "TTreeIndex.h"
17 
18 #include "TTreeFormula.h"
19 #include "TTree.h"
20 #include "TBuffer.h"
21 #include "TMath.h"
22 
24 
25 
26 struct IndexSortComparator {
27 
28  IndexSortComparator(Long64_t *major, Long64_t *minor)
29  : fValMajor(major), fValMinor(minor)
30  {}
31 
32  template<typename Index>
33  bool operator()(Index i1, Index i2) {
34  if( *(fValMajor + i1) == *(fValMajor + i2) )
35  return *(fValMinor + i1) < *(fValMinor + i2);
36  else
37  return *(fValMajor + i1) < *(fValMajor + i2);
38  }
39 
40  // pointers to the start of index values tables keeping upper 64bit and lower 64bit
41  // of combined indexed 128bit value
42  Long64_t *fValMajor, *fValMinor;
43 };
44 
45 
46 ////////////////////////////////////////////////////////////////////////////////
47 /// Default constructor for TTreeIndex
48 
50 {
51  fTree = 0;
52  fN = 0;
53  fIndexValues = 0;
55  fIndex = 0;
56  fMajorFormula = 0;
57  fMinorFormula = 0;
60 }
61 
62 ////////////////////////////////////////////////////////////////////////////////
63 /// Normal constructor for TTreeIndex
64 ///
65 /// Build an index table using the leaves of Tree T with major & minor names
66 /// The index is built with the expressions given in "majorname" and "minorname".
67 ///
68 /// a Long64_t array fIndexValues is built with:
69 ///
70 /// - major = the value of majorname converted to an integer
71 /// - minor = the value of minorname converted to an integer
72 /// - fIndexValues[i] = major<<31 + minor
73 ///
74 /// This array is sorted. The sorted fIndex[i] contains the serial number
75 /// in the Tree corresponding to the pair "major,minor" in fIndexvalues[i].
76 ///
77 /// Once the index is computed, one can retrieve one entry via
78 /// ~~~{.cpp}
79 /// T->GetEntryWithIndex(majornumber, minornumber)
80 /// ~~~
81 /// Example:
82 /// ~~~{.cpp}
83 /// tree.BuildIndex("Run","Event"); //creates an index using leaves Run and Event
84 /// tree.GetEntryWithIndex(1234,56789); // reads entry corresponding to
85 /// // Run=1234 and Event=56789
86 /// ~~~
87 /// Note that majorname and minorname may be expressions using original
88 /// Tree variables eg: "run-90000", "event +3*xx". However the result
89 /// must be integer.
90 ///
91 /// In case an expression is specified, the equivalent expression must be computed
92 /// when calling GetEntryWithIndex.
93 ///
94 /// To build an index with only majorname, specify minorname="0" (default)
95 ///
96 /// ## TreeIndex and Friend Trees
97 ///
98 /// Assuming a parent Tree T and a friend Tree TF, the following cases are supported:
99 /// - CASE 1: T->GetEntry(entry) is called
100 /// In this case, the serial number entry is used to retrieve
101 /// the data in both Trees.
102 /// - CASE 2: T->GetEntry(entry) is called, TF has a TreeIndex
103 /// the expressions given in major/minorname of TF are used
104 /// to compute the value pair major,minor with the data in T.
105 /// TF->GetEntryWithIndex(major,minor) is then called (tricky case!)
106 /// - CASE 3: T->GetEntryWithIndex(major,minor) is called.
107 /// It is assumed that both T and TF have a TreeIndex built using
108 /// the same major and minor name.
109 ///
110 /// ## Saving the TreeIndex
111 ///
112 /// Once the index is built, it can be saved with the TTree object
113 /// with tree.Write(); (if the file has been open in "update" mode).
114 ///
115 /// The most convenient place to create the index is at the end of
116 /// the filling process just before saving the Tree header.
117 /// If a previous index was computed, it is redefined by this new call.
118 ///
119 /// Note that this function can also be applied to a TChain.
120 ///
121 /// The return value is the number of entries in the Index (< 0 indicates failure)
122 ///
123 /// It is possible to play with different TreeIndex in the same Tree.
124 /// see comments in TTree::SetTreeIndex.
125 
126 TTreeIndex::TTreeIndex(const TTree *T, const char *majorname, const char *minorname)
127  : TVirtualIndex()
128 {
129  fTree = (TTree*)T;
130  fN = 0;
131  fIndexValues = 0;
132  fIndexValuesMinor = 0;
133  fIndex = 0;
134  fMajorFormula = 0;
135  fMinorFormula = 0;
138  fMajorName = majorname;
139  fMinorName = minorname;
140  if (!T) return;
141  fN = T->GetEntries();
142  if (fN <= 0) {
143  MakeZombie();
144  Error("TreeIndex","Cannot build a TreeIndex with a Tree having no entries");
145  return;
146  }
147 
148  GetMajorFormula();
149  GetMinorFormula();
150  if (!fMajorFormula || !fMinorFormula) {
151  MakeZombie();
152  Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
153  return;
154  }
155  if ((fMajorFormula->GetNdim() != 1) || (fMinorFormula->GetNdim() != 1)) {
156  MakeZombie();
157  Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
158  return;
159  }
160  // accessing array elements should be OK
161  //if ((fMajorFormula->GetMultiplicity() != 0) || (fMinorFormula->GetMultiplicity() != 0)) {
162  // MakeZombie();
163  // Error("TreeIndex","Cannot build the index with major=%s, minor=%s that cannot be arrays",fMajorName.Data(), fMinorName.Data());
164  // return;
165  //}
166 
167  Long64_t *tmp_major = new Long64_t[fN];
168  Long64_t *tmp_minor = new Long64_t[fN];
169  Long64_t i;
170  Long64_t oldEntry = fTree->GetReadEntry();
171  Int_t current = -1;
172  for (i=0;i<fN;i++) {
173  Long64_t centry = fTree->LoadTree(i);
174  if (centry < 0) break;
175  if (fTree->GetTreeNumber() != current) {
176  current = fTree->GetTreeNumber();
179  }
180  tmp_major[i] = (Long64_t) fMajorFormula->EvalInstance<LongDouble_t>();
181  tmp_minor[i] = (Long64_t) fMinorFormula->EvalInstance<LongDouble_t>();
182  }
183  fIndex = new Long64_t[fN];
184  for(i = 0; i < fN; i++) { fIndex[i] = i; }
185  std::sort(fIndex, fIndex + fN, IndexSortComparator(tmp_major, tmp_minor) );
186  //TMath::Sort(fN,w,fIndex,0);
187  fIndexValues = new Long64_t[fN];
189  for (i=0;i<fN;i++) {
190  fIndexValues[i] = tmp_major[fIndex[i]];
191  fIndexValuesMinor[i] = tmp_minor[fIndex[i]];
192  }
193 
194  delete [] tmp_major;
195  delete [] tmp_minor;
196  fTree->LoadTree(oldEntry);
197 }
198 
199 ////////////////////////////////////////////////////////////////////////////////
200 /// Destructor.
201 
203 {
204  if (fTree && fTree->GetTreeIndex() == this) fTree->SetTreeIndex(0);
205  delete [] fIndexValues; fIndexValues = 0;
206  delete [] fIndexValuesMinor; fIndexValuesMinor = 0;
207  delete [] fIndex; fIndex = 0;
208  delete fMajorFormula; fMajorFormula = 0;
209  delete fMinorFormula; fMinorFormula = 0;
212 }
213 
214 ////////////////////////////////////////////////////////////////////////////////
215 /// Append 'add' to this index. Entry 0 in add will become entry n+1 in this.
216 /// If delaySort is true, do not sort the value, then you must call
217 /// Append(0,kFALSE);
218 
219 void TTreeIndex::Append(const TVirtualIndex *add, Bool_t delaySort )
220 {
221 
222  if (add && add->GetN()) {
223  // Create new buffer (if needed)
224 
225  const TTreeIndex *ti_add = dynamic_cast<const TTreeIndex*>(add);
226  if (ti_add == 0) {
227  Error("Append","Can only Append a TTreeIndex to a TTreeIndex but got a %s",
228  add->IsA()->GetName());
229  }
230 
231  Long64_t oldn = fN;
232  fN += add->GetN();
233 
234  Long64_t *oldIndex = fIndex;
235  Long64_t *oldValues = GetIndexValues();
236  Long64_t *oldValues2 = GetIndexValuesMinor();
237 
238  fIndex = new Long64_t[fN];
239  fIndexValues = new Long64_t[fN];
241 
242  // Copy data
243  Long_t size = sizeof(Long64_t) * oldn;
244  Long_t add_size = sizeof(Long64_t) * add->GetN();
245 
246  memcpy(fIndex,oldIndex, size);
247  memcpy(fIndexValues,oldValues, size);
248  memcpy(fIndexValuesMinor,oldValues2, size);
249 
250  Long64_t *addIndex = ti_add->GetIndex();
251  Long64_t *addValues = ti_add->GetIndexValues();
252  Long64_t *addValues2 = ti_add->GetIndexValuesMinor();
253 
254  memcpy(fIndex + oldn, addIndex, add_size);
255  memcpy(fIndexValues + oldn, addValues, add_size);
256  memcpy(fIndexValuesMinor + oldn, addValues2, add_size);
257  for(Int_t i = 0; i < add->GetN(); i++) {
258  fIndex[oldn + i] += oldn;
259  }
260 
261  delete [] oldIndex;
262  delete [] oldValues;
263  delete [] oldValues2;
264  }
265 
266  // Sort.
267  if (!delaySort) {
268  Long64_t *addValues = GetIndexValues();
269  Long64_t *addValues2 = GetIndexValuesMinor();
270  Long64_t *ind = fIndex;
271  Long64_t *conv = new Long64_t[fN];
272 
273  for(Long64_t i = 0; i < fN; i++) { conv[i] = i; }
274  std::sort(conv, conv+fN, IndexSortComparator(addValues, addValues2) );
275  //Long64_t *w = fIndexValues;
276  //TMath::Sort(fN,w,conv,0);
277 
278  fIndex = new Long64_t[fN];
279  fIndexValues = new Long64_t[fN];
281 
282  for (Int_t i=0;i<fN;i++) {
283  fIndex[i] = ind[conv[i]];
284  fIndexValues[i] = addValues[conv[i]];
285  fIndexValuesMinor[i] = addValues2[conv[i]];
286  }
287  delete [] addValues;
288  delete [] addValues2;
289  delete [] ind;
290  delete [] conv;
291  }
292 }
293 
294 
295 
296 ////////////////////////////////////////////////////////////////////////////////
297 /// conversion from old 64bit indexes
298 /// return true if index was converted
299 
301 {
302  if( !fIndexValuesMinor && fN ) {
304  for(int i=0; i<fN; i++) {
305  fIndexValuesMinor[i] = (fIndexValues[i] & 0x7fffffff);
306  fIndexValues[i] >>= 31;
307  }
308  return true;
309  }
310  return false;
311 }
312 
313 
314 
315 ////////////////////////////////////////////////////////////////////////////////
316 /// Returns the entry number in this (friend) Tree corresponding to entry in
317 /// the master Tree 'parent'.
318 /// In case this (friend) Tree and 'master' do not share an index with the same
319 /// major and minor name, the entry serial number in the (friend) tree
320 /// and in the master Tree are assumed to be the same
321 
323 {
324  if (!parent) return -3;
325  // We reached the end of the parent tree
326  Long64_t pentry = parent->GetReadEntry();
327  if (pentry >= parent->GetEntries())
328  return -2;
329  GetMajorFormulaParent(parent);
330  GetMinorFormulaParent(parent);
331  if (!fMajorFormulaParent || !fMinorFormulaParent) return -1;
333  // The Tree Index in the friend has a pair majorname,minorname
334  // not available in the parent Tree T.
335  // if the friend Tree has less entries than the parent, this is an error
336  if (pentry >= fTree->GetEntries()) return -2;
337  // otherwise we ignore the Tree Index and return the entry number
338  // in the parent Tree.
339  return pentry;
340  }
341 
342  // majorname, minorname exist in the parent Tree
343  // we find the current values pair majorv,minorv in the parent Tree
346  Long64_t majorv = (Long64_t)majord;
347  Long64_t minorv = (Long64_t)minord;
348  // we check if this pair exist in the index.
349  // if yes, we return the corresponding entry number
350  // if not the function returns -1
351  return fTree->GetEntryNumberWithIndex(majorv,minorv);
352 }
353 
354 
355 ////////////////////////////////////////////////////////////////////////////////
356 /// find position where major|minor values are in the IndexValues tables
357 /// this is the index in IndexValues table, not entry# !
358 /// use lower_bound STD algorithm.
359 
361 {
362  Long64_t mid, step, pos = 0, count = fN;
363  // find lower bound using bisection
364  while( count > 0 ) {
365  step = count / 2;
366  mid = pos + step;
367  // check if *mid < major|minor
368  if( fIndexValues[mid] < major
369  || ( fIndexValues[mid] == major && fIndexValuesMinor[mid] < minor ) ) {
370  pos = mid+1;
371  count -= step + 1;
372  } else
373  count = step;
374  }
375  return pos;
376 }
377 
378 
379 ////////////////////////////////////////////////////////////////////////////////
380 /// Return entry number corresponding to major and minor number.
381 /// Note that this function returns only the entry number, not the data
382 /// To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
383 /// the BuildIndex function has created a table of Double_t* of sorted values
384 /// corresponding to val = major<<31 + minor;
385 /// The function performs binary search in this sorted table.
386 /// If it finds a pair that maches val, it returns directly the
387 /// index in the table.
388 /// If an entry corresponding to major and minor is not found, the function
389 /// returns the index of the major,minor pair immediately lower than the
390 /// requested value, ie it will return -1 if the pair is lower than
391 /// the first entry in the index.
392 ///
393 /// See also GetEntryNumberWithIndex
394 
396 {
397  if (fN == 0) return -1;
398 
399  Long64_t pos = FindValues(major, minor);
400  if( pos < fN && fIndexValues[pos] == major && fIndexValuesMinor[pos] == minor )
401  return fIndex[pos];
402  if( --pos < 0 )
403  return -1;
404  return fIndex[pos];
405 }
406 
407 
408 ////////////////////////////////////////////////////////////////////////////////
409 /// Return entry number corresponding to major and minor number.
410 /// Note that this function returns only the entry number, not the data
411 /// To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
412 /// the BuildIndex function has created a table of Double_t* of sorted values
413 /// corresponding to val = major<<31 + minor;
414 /// The function performs binary search in this sorted table.
415 /// If it finds a pair that maches val, it returns directly the
416 /// index in the table, otherwise it returns -1.
417 ///
418 /// See also GetEntryNumberWithBestIndex
419 
421 {
422  if (fN == 0) return -1;
423 
424  Long64_t pos = FindValues(major, minor);
425  if( pos < fN && fIndexValues[pos] == major && fIndexValuesMinor[pos] == minor )
426  return fIndex[pos];
427  return -1;
428 }
429 
430 
431 ////////////////////////////////////////////////////////////////////////////////
432 
434 {
435  return fIndexValuesMinor;
436 }
437 
438 
439 
440 ////////////////////////////////////////////////////////////////////////////////
441 /// Return a pointer to the TreeFormula corresponding to the majorname.
442 
444 {
445  if (!fMajorFormula) {
448  }
449  return fMajorFormula;
450 }
451 
452 ////////////////////////////////////////////////////////////////////////////////
453 /// Return a pointer to the TreeFormula corresponding to the minorname.
454 
456 {
457  if (!fMinorFormula) {
460  }
461  return fMinorFormula;
462 }
463 
464 ////////////////////////////////////////////////////////////////////////////////
465 /// Return a pointer to the TreeFormula corresponding to the majorname in parent tree.
466 
468 {
469  if (!fMajorFormulaParent) {
470  // Prevent TTreeFormula from finding any of the branches in our TTree even if it
471  // is a friend of the parent TTree.
473  fMajorFormulaParent = new TTreeFormula("MajorP",fMajorName.Data(),const_cast<TTree*>(parent));
475  }
476  if (fMajorFormulaParent->GetTree() != parent) {
477  fMajorFormulaParent->SetTree(const_cast<TTree*>(parent));
479  }
480  return fMajorFormulaParent;
481 }
482 
483 ////////////////////////////////////////////////////////////////////////////////
484 /// Return a pointer to the TreeFormula corresponding to the minorname in parent tree.
485 
487 {
488  if (!fMinorFormulaParent) {
489  // Prevent TTreeFormula from finding any of the branches in our TTree even if it
490  // is a friend of the parent TTree.
492  fMinorFormulaParent = new TTreeFormula("MinorP",fMinorName.Data(),const_cast<TTree*>(parent));
494  }
495  if (fMinorFormulaParent->GetTree() != parent) {
496  fMinorFormulaParent->SetTree(const_cast<TTree*>(parent));
498  }
499  return fMinorFormulaParent;
500 }
501 
502 ////////////////////////////////////////////////////////////////////////////////
503 /// Return kTRUE if index can be applied to the TTree
504 
506 {
507  auto *majorFormula = GetMajorFormulaParent(parent);
508  auto *minorFormula = GetMinorFormulaParent(parent);
509  if ((majorFormula == nullptr || majorFormula->GetNdim() == 0) ||
510  (minorFormula == nullptr || minorFormula->GetNdim() == 0))
511  return kFALSE;
512  return kTRUE;
513 }
514 
515 ////////////////////////////////////////////////////////////////////////////////
516 /// Print the table with : serial number, majorname, minorname.
517 /// - if option = "10" print only the first 10 entries
518 /// - if option = "100" print only the first 100 entries
519 /// - if option = "1000" print only the first 1000 entries
520 
521 void TTreeIndex::Print(Option_t * option) const
522 {
523  TString opt = option;
524  Bool_t printEntry = kFALSE;
525  Long64_t n = fN;
526  if (opt.Contains("10")) n = 10;
527  if (opt.Contains("100")) n = 100;
528  if (opt.Contains("1000")) n = 1000;
529  if (opt.Contains("all")) {
530  printEntry = kTRUE;
531  }
532 
533  if (printEntry) {
534  Printf("\n*****************************************************************");
535  Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
536  Printf("*****************************************************************");
537  Printf("%8s : %16s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data(),"entry number");
538  Printf("*****************************************************************");
539  for (Long64_t i=0;i<n;i++) {
540  Printf("%8lld : %8lld : %8lld : %8lld",
541  i, fIndexValues[i], GetIndexValuesMinor()[i], fIndex[i]);
542  }
543 
544  } else {
545  Printf("\n**********************************************");
546  Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
547  Printf("**********************************************");
548  Printf("%8s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data());
549  Printf("**********************************************");
550  for (Long64_t i=0;i<n;i++) {
551  Printf("%8lld : %8lld : %8lld",
552  i, fIndexValues[i],GetIndexValuesMinor()[i]);
553  }
554  }
555 }
556 
557 ////////////////////////////////////////////////////////////////////////////////
558 /// Stream an object of class TTreeIndex.
559 /// Note that this Streamer should be changed to an automatic Streamer
560 /// once TStreamerInfo supports an index of type Long64_t
561 
562 void TTreeIndex::Streamer(TBuffer &R__b)
563 {
564  UInt_t R__s, R__c;
565  if (R__b.IsReading()) {
566  Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
567  TVirtualIndex::Streamer(R__b);
568  fMajorName.Streamer(R__b);
569  fMinorName.Streamer(R__b);
570  R__b >> fN;
571  fIndexValues = new Long64_t[fN];
573  if( R__v > 1 ) {
576  } else {
577  ConvertOldToNew();
578  }
579  fIndex = new Long64_t[fN];
580  R__b.ReadFastArray(fIndex,fN);
581  R__b.CheckByteCount(R__s, R__c, TTreeIndex::IsA());
582  } else {
583  R__c = R__b.WriteVersion(TTreeIndex::IsA(), kTRUE);
584  TVirtualIndex::Streamer(R__b);
585  fMajorName.Streamer(R__b);
586  fMinorName.Streamer(R__b);
587  R__b << fN;
590  R__b.WriteFastArray(fIndex, fN);
591  R__b.SetByteCount(R__c, kTRUE);
592  }
593 }
594 
595 ////////////////////////////////////////////////////////////////////////////////
596 /// Called by TChain::LoadTree when the parent chain changes it's tree.
597 
599 {
602  if (fMajorFormulaParent) {
603  if (parent) fMajorFormulaParent->SetTree(const_cast<TTree*>(parent));
605  }
606  if (fMinorFormulaParent) {
607  if (parent) fMinorFormulaParent->SetTree(const_cast<TTree*>(parent));
609  }
610 }
611 ////////////////////////////////////////////////////////////////////////////////
612 /// this function is called by TChain::LoadTree and TTreePlayer::UpdateFormulaLeaves
613 /// when a new Tree is loaded.
614 /// Because Trees in a TChain may have a different list of leaves, one
615 /// must update the leaves numbers in the TTreeFormula used by the TreeIndex.
616 
618 {
619  fTree = (TTree*)T;
620 }
621 
n
const Int_t n
Definition: legend1.C:16
TTreeFormula::SetQuickLoad
void SetQuickLoad(Bool_t quick)
Definition: TTreeFormula.h:207
TTreeIndex::Print
virtual void Print(Option_t *option="") const
Print the table with : serial number, majorname, minorname.
Definition: TTreeIndex.cxx:521
kTRUE
const Bool_t kTRUE
Definition: RtypesCore.h:100
TTreeIndex.h
Version_t
short Version_t
Definition: RtypesCore.h:65
Option_t
const char Option_t
Definition: RtypesCore.h:66
TBuffer::WriteVersion
virtual UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE)=0
TBuffer::SetByteCount
virtual void SetByteCount(UInt_t cntpos, Bool_t packInVersion=kFALSE)=0
TTreeIndex::FindValues
Long64_t FindValues(Long64_t major, Long64_t minor) const
find position where major|minor values are in the IndexValues tables this is the index in IndexValues...
Definition: TTreeIndex.cxx:360
TVirtualIndex
Abstract interface for Tree Index.
Definition: TVirtualIndex.h:30
TString::Data
const char * Data() const
Definition: TString.h:369
TBuffer::ReadFastArray
virtual void ReadFastArray(Bool_t *b, Int_t n)=0
ClassImp
#define ClassImp(name)
Definition: Rtypes.h:364
TNamed::GetTitle
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
TTreeIndex::GetEntryNumberWithBestIndex
virtual Long64_t GetEntryNumberWithBestIndex(Long64_t major, Long64_t minor) const
Return entry number corresponding to major and minor number.
Definition: TTreeIndex.cxx:395
TTree::kFindBranch
@ kFindBranch
Definition: TTree.h:208
TTreeIndex::GetMajorFormulaParent
TTreeFormula * GetMajorFormulaParent(const TTree *parent)
Return a pointer to the TreeFormula corresponding to the majorname in parent tree.
Definition: TTreeIndex.cxx:467
Long64_t
long long Long64_t
Definition: RtypesCore.h:80
TObject::Error
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:893
TTree
A TTree represents a columnar dataset.
Definition: TTree.h:79
TTreeFormula.h
TTree::SetTreeIndex
virtual void SetTreeIndex(TVirtualIndex *index)
The current TreeIndex is replaced by the new index.
Definition: TTree.cxx:9201
TTreeIndex::GetIndexValues
virtual Long64_t * GetIndexValues() const
Definition: TTreeIndex.h:61
TTree::kFindLeaf
@ kFindLeaf
Definition: TTree.h:209
TTreeIndex::TTreeIndex
TTreeIndex()
Default constructor for TTreeIndex.
Definition: TTreeIndex.cxx:49
TString::Contains
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:624
operator()
TRObject operator()(const T1 &t1) const
Definition: TRFunctionImport__oprtr.h:14
TTreeIndex::~TTreeIndex
virtual ~TTreeIndex()
Destructor.
Definition: TTreeIndex.cxx:202
TBuffer
Buffer base class used for serializing objects.
Definition: TBuffer.h:43
TTree.h
TBuffer::CheckByteCount
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
TString
Basic string class.
Definition: TString.h:136
TTreeIndex::SetTree
virtual void SetTree(const TTree *T)
this function is called by TChain::LoadTree and TTreePlayer::UpdateFormulaLeaves when a new Tree is l...
Definition: TTreeIndex.cxx:617
TTreeIndex::GetMinorFormula
virtual TTreeFormula * GetMinorFormula()
Return a pointer to the TreeFormula corresponding to the minorname.
Definition: TTreeIndex.cxx:455
TTreeIndex::Append
virtual void Append(const TVirtualIndex *, Bool_t delaySort=kFALSE)
Append 'add' to this index.
Definition: TTreeIndex.cxx:219
TTreeIndex::fMajorFormula
TTreeFormula * fMajorFormula
! Pointer to major TreeFormula
Definition: TTreeIndex.h:38
bool
TTreeIndex::ConvertOldToNew
bool ConvertOldToNew()
conversion from old 64bit indexes return true if index was converted
Definition: TTreeIndex.cxx:300
TTreeFormula::UpdateFormulaLeaves
virtual void UpdateFormulaLeaves()
This function is called TTreePlayer::UpdateFormulaLeaves, itself called by TChain::LoadTree when a ne...
Definition: TTreeFormula.cxx:5066
TTreeFormula::SetTree
virtual void SetTree(TTree *tree)
Definition: TTreeFormula.h:208
TBuffer.h
TTreeIndex::fIndexValues
Long64_t * fIndexValues
[fN] Sorted index values, higher 64bits
Definition: TTreeIndex.h:35
size
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
TBuffer::WriteFastArray
virtual void WriteFastArray(const Bool_t *b, Int_t n)=0
TTree::GetReadEntry
virtual Long64_t GetReadEntry() const
Definition: TTree.h:504
TTree::LoadTree
virtual Long64_t LoadTree(Long64_t entry)
Set current entry.
Definition: TTree.cxx:6410
TTreeIndex::fIndexValuesMinor
Long64_t * fIndexValuesMinor
[fN] Sorted index values, lower 64bits
Definition: TTreeIndex.h:36
kFALSE
const Bool_t kFALSE
Definition: RtypesCore.h:101
TTreeIndex::GetEntryNumberWithIndex
virtual Long64_t GetEntryNumberWithIndex(Long64_t major, Long64_t minor) const
Return entry number corresponding to major and minor number.
Definition: TTreeIndex.cxx:420
TTreeFormula::EvalInstance
T EvalInstance(Int_t i=0, const char *stringStack[]=0)
Evaluate this treeformula.
Definition: TTreeFormula.cxx:3936
Long_t
long Long_t
Definition: RtypesCore.h:54
TTreeFormula
Used to pass a selection expression to the Tree drawing routine.
Definition: TTreeFormula.h:58
TTreeIndex::fMajorName
TString fMajorName
Index major name.
Definition: TTreeIndex.h:32
TTreeIndex::UpdateFormulaLeaves
virtual void UpdateFormulaLeaves(const TTree *parent)
Called by TChain::LoadTree when the parent chain changes it's tree.
Definition: TTreeIndex.cxx:598
UInt_t
unsigned int UInt_t
Definition: RtypesCore.h:46
TObject::MakeZombie
void MakeZombie()
Definition: TObject.h:49
TTreeIndex::fMinorFormula
TTreeFormula * fMinorFormula
! Pointer to minor TreeFormula
Definition: TTreeIndex.h:39
TTree::GetTreeNumber
virtual Int_t GetTreeNumber() const
Definition: TTree.h:514
TTree::kGetBranch
@ kGetBranch
Definition: TTree.h:211
TTree::TFriendLock
Helper class to prevent infinite recursion in the usage of TTree Friends.
Definition: TTree.h:184
TBuffer::ReadVersion
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
TTreeIndex::fMinorName
TString fMinorName
Index minor name.
Definition: TTreeIndex.h:33
Printf
void Printf(const char *fmt,...)
TTree::GetEntryNumberWithIndex
virtual Long64_t GetEntryNumberWithIndex(Long64_t major, Long64_t minor=0) const
Return entry number corresponding to major and minor number.
Definition: TTree.cxx:5847
TTreeIndex
A Tree Index with majorname and minorname.
Definition: TTreeIndex.h:29
TBuffer::IsReading
Bool_t IsReading() const
Definition: TBuffer.h:86
TTree::GetTreeIndex
virtual TVirtualIndex * GetTreeIndex() const
Definition: TTree.h:513
TTreeIndex::IsValidFor
virtual Bool_t IsValidFor(const TTree *parent)
Return kTRUE if index can be applied to the TTree.
Definition: TTreeIndex.cxx:505
TTreeIndex::GetMinorFormulaParent
TTreeFormula * GetMinorFormulaParent(const TTree *parent)
Return a pointer to the TreeFormula corresponding to the minorname in parent tree.
Definition: TTreeIndex.cxx:486
LongDouble_t
long double LongDouble_t
Definition: RtypesCore.h:61
Double_t
double Double_t
Definition: RtypesCore.h:59
TTree::kGetLeaf
@ kGetLeaf
Definition: TTree.h:216
TTreeIndex::fMajorFormulaParent
TTreeFormula * fMajorFormulaParent
! Pointer to major TreeFormula in Parent tree (if any)
Definition: TTreeIndex.h:40
ROOT::Math::Chebyshev::T
double T(double x)
Definition: ChebyshevPol.h:34
ROOT::v5::TFormula::GetNdim
virtual Int_t GetNdim() const
Definition: TFormula.h:237
TNamed::GetName
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
RooFit::Index
RooCmdArg Index(RooCategory &icat)
Definition: RooGlobalFunc.cxx:107
TTreeIndex::GetIndex
virtual Long64_t * GetIndex() const
Definition: TTreeIndex.h:60
TTreeIndex::fIndex
Long64_t * fIndex
[fN] Index of sorted values
Definition: TTreeIndex.h:37
TTreeIndex::GetEntryNumberFriend
virtual Long64_t GetEntryNumberFriend(const TTree *parent)
Returns the entry number in this (friend) Tree corresponding to entry in the master Tree 'parent'.
Definition: TTreeIndex.cxx:322
TTree::GetEntries
virtual Long64_t GetEntries() const
Definition: TTree.h:458
TTreeIndex::fN
Long64_t fN
Number of entries.
Definition: TTreeIndex.h:34
TTreeIndex::GetIndexValuesMinor
virtual Long64_t * GetIndexValuesMinor() const
Definition: TTreeIndex.cxx:433
TTreeIndex::fMinorFormulaParent
TTreeFormula * fMinorFormulaParent
! Pointer to minor TreeFormula in Parent tree (if any)
Definition: TTreeIndex.h:41
TTreeFormula::GetTree
virtual TTree * GetTree() const
Definition: TTreeFormula.h:210
TVirtualIndex::fTree
TTree * fTree
Definition: TVirtualIndex.h:33
TMath.h
TVirtualIndex::GetN
virtual Long64_t GetN() const =0
int
TTreeIndex::GetMajorFormula
virtual TTreeFormula * GetMajorFormula()
Return a pointer to the TreeFormula corresponding to the majorname.
Definition: TTreeIndex.cxx:443