Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
13A 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
23#include <cstring> // std::strlen
24
26
27
29
31 : fValMajor(major), fValMinor(minor)
32 {}
33
34 template<typename Index>
35 bool operator()(Index i1, Index i2) {
36 if( *(fValMajor + i1) == *(fValMajor + i2) )
37 return *(fValMinor + i1) < *(fValMinor + i2);
38 else
39 return *(fValMajor + i1) < *(fValMajor + i2);
40 }
41
42 // pointers to the start of index values tables keeping upper 64bit and lower 64bit
43 // of combined indexed 128bit value
45};
46
47
48////////////////////////////////////////////////////////////////////////////////
49/// Default constructor for TTreeIndex
50
52{
53 fTree = nullptr;
54 fN = 0;
55 fIndexValues = nullptr;
56 fIndexValuesMinor = nullptr;
57 fIndex = nullptr;
58 fMajorFormula = nullptr;
59 fMinorFormula = nullptr;
60 fMajorFormulaParent = nullptr;
61 fMinorFormulaParent = nullptr;
62}
63
64////////////////////////////////////////////////////////////////////////////////
65/// Normal constructor for TTreeIndex
66///
67/// Build an index table using the leaves of Tree T with major & minor names
68/// The index is built with the expressions given in "majorname" and "minorname".
69///
70/// a Long64_t array fIndexValues is built with:
71///
72/// - major = the value of majorname converted to an integer
73/// - minor = the value of minorname converted to an integer
74/// - fIndexValues[i] = major<<31 + minor
75///
76/// This array is sorted. The sorted fIndex[i] contains the serial number
77/// in the Tree corresponding to the pair "major,minor" in fIndexvalues[i].
78///
79/// Once the index is computed, one can retrieve one entry via
80/// ~~~{.cpp}
81/// T->GetEntryWithIndex(majornumber, minornumber)
82/// ~~~
83/// Example:
84/// ~~~{.cpp}
85/// tree.BuildIndex("Run","Event"); //creates an index using leaves Run and Event
86/// tree.GetEntryWithIndex(1234,56789); // reads entry corresponding to
87/// // Run=1234 and Event=56789
88/// ~~~
89/// Note that majorname and minorname may be expressions using original
90/// Tree variables eg: "run-90000", "event +3*xx". However the result
91/// must be integer.
92///
93/// In case an expression is specified, the equivalent expression must be computed
94/// when calling GetEntryWithIndex.
95///
96/// To build an index with only majorname, specify minorname="0" (default)
97///
98/// ## TreeIndex and Friend Trees
99///
100/// Assuming a parent Tree T and a friend Tree TF, the following cases are supported:
101/// - CASE 1: T->GetEntry(entry) is called
102/// In this case, the serial number entry is used to retrieve
103/// the data in both Trees.
104/// - CASE 2: T->GetEntry(entry) is called, TF has a TreeIndex
105/// the expressions given in major/minorname of TF are used
106/// to compute the value pair major,minor with the data in T.
107/// TF->GetEntryWithIndex(major,minor) is then called (tricky case!)
108/// - CASE 3: T->GetEntryWithIndex(major,minor) is called.
109/// It is assumed that both T and TF have a TreeIndex built using
110/// the same major and minor name.
111///
112/// ## Saving the TreeIndex
113///
114/// Once the index is built, it can be saved with the TTree object
115/// with tree.Write(); (if the file has been open in "update" mode).
116///
117/// The most convenient place to create the index is at the end of
118/// the filling process just before saving the Tree header.
119/// If a previous index was computed, it is redefined by this new call.
120///
121/// Note that this function can also be applied to a TChain.
122///
123/// The return value is the number of entries in the Index (< 0 indicates failure)
124///
125/// It is possible to play with different TreeIndex in the same Tree.
126/// see comments in TTree::SetTreeIndex.
127
128TTreeIndex::TTreeIndex(const TTree *T, const char *majorname, const char *minorname)
129 : TVirtualIndex()
130{
131 fTree = (TTree*)T;
132 fN = 0;
133 fIndexValues = nullptr;
134 fIndexValuesMinor = nullptr;
135 fIndex = nullptr;
136 fMajorFormula = nullptr;
137 fMinorFormula = nullptr;
138 fMajorFormulaParent = nullptr;
139 fMinorFormulaParent = nullptr;
140 fMajorName = majorname;
141 fMinorName = minorname;
142 if (!T) return;
143 fN = T->GetEntries();
144 if (fN <= 0) {
145 MakeZombie();
146 Error("TreeIndex","Cannot build a TreeIndex with a Tree having no entries");
147 return;
148 }
149
152 if (!fMajorFormula || !fMinorFormula) {
153 MakeZombie();
154 Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
155 return;
156 }
157 if ((fMajorFormula->GetNdim() != 1) || (fMinorFormula->GetNdim() != 1)) {
158 MakeZombie();
159 Error("TreeIndex","Cannot build the index with major=%s, minor=%s",fMajorName.Data(), fMinorName.Data());
160 return;
161 }
162 // accessing array elements should be OK
163 //if ((fMajorFormula->GetMultiplicity() != 0) || (fMinorFormula->GetMultiplicity() != 0)) {
164 // MakeZombie();
165 // Error("TreeIndex","Cannot build the index with major=%s, minor=%s that cannot be arrays",fMajorName.Data(), fMinorName.Data());
166 // return;
167 //}
168
169 Long64_t *tmp_major = new Long64_t[fN];
170 Long64_t *tmp_minor = new Long64_t[fN];
171 Long64_t i;
172 Long64_t oldEntry = fTree->GetReadEntry();
173 Int_t current = -1;
174 for (i=0;i<fN;i++) {
175 Long64_t centry = fTree->LoadTree(i);
176 if (centry < 0) break;
177 if (fTree->GetTreeNumber() != current) {
178 current = fTree->GetTreeNumber();
181 }
182 auto GetAndRangeCheck = [this](bool isMajor, Long64_t entry) {
183 LongDouble_t ret = (isMajor ? fMajorFormula : fMinorFormula)->EvalInstance<LongDouble_t>();
184 // Check whether the value (vs significant bits) of ldRet can represent
185 // the full precision of the returned value. If we return 10^60, the
186 // value fits into a long double, but if sizeof(long double) ==
187 // sizeof(double) it cannot store the ones: the value returned by
188 // EvalInstance() only stores the higher bits.
189 LongDouble_t retCloserToZero = ret;
190 if (ret > 0)
191 retCloserToZero -= 1;
192 else
193 retCloserToZero += 1;
194 if (retCloserToZero == ret) {
195 Warning("TTreeIndex",
196 "In tree entry %lld, %s value %s=%Lf possibly out of range for internal `long double`", entry,
197 isMajor ? "major" : "minor", isMajor ? fMajorName.Data() : fMinorName.Data(), ret);
198 }
199 return ret;
200 };
201 tmp_major[i] = GetAndRangeCheck(true, i);
202 tmp_minor[i] = GetAndRangeCheck(false, i);
203 }
204 fIndex = new Long64_t[fN];
205 for(i = 0; i < fN; i++) { fIndex[i] = i; }
206 std::sort(fIndex, fIndex + fN, IndexSortComparator(tmp_major, tmp_minor) );
207 //TMath::Sort(fN,w,fIndex,0);
208 fIndexValues = new Long64_t[fN];
210 for (i=0;i<fN;i++) {
211 fIndexValues[i] = tmp_major[fIndex[i]];
212 fIndexValuesMinor[i] = tmp_minor[fIndex[i]];
213 }
214
215 delete [] tmp_major;
216 delete [] tmp_minor;
217 fTree->LoadTree(oldEntry);
218}
219
220////////////////////////////////////////////////////////////////////////////////
221/// Destructor.
222
224{
225 if (fTree && fTree->GetTreeIndex() == this) fTree->SetTreeIndex(nullptr);
226 delete [] fIndexValues; fIndexValues = nullptr;
227 delete [] fIndexValuesMinor; fIndexValuesMinor = nullptr;
228 delete [] fIndex; fIndex = nullptr;
229 delete fMajorFormula; fMajorFormula = nullptr;
230 delete fMinorFormula; fMinorFormula = nullptr;
233}
234
235////////////////////////////////////////////////////////////////////////////////
236/// Append 'add' to this index. Entry 0 in add will become entry n+1 in this.
237/// If delaySort is true, do not sort the value, then you must call
238/// Append(0,false);
239
240void TTreeIndex::Append(const TVirtualIndex *add, bool delaySort )
241{
242
243 if (add && add->GetN()) {
244 // Create new buffer (if needed)
245
246 const TTreeIndex *ti_add = dynamic_cast<const TTreeIndex*>(add);
247 if (ti_add == nullptr) {
248 Error("Append","Can only Append a TTreeIndex to a TTreeIndex but got a %s",
249 add->IsA()->GetName());
250 }
251
252 Long64_t oldn = fN;
253 fN += add->GetN();
254
255 Long64_t *oldIndex = fIndex;
256 Long64_t *oldValues = GetIndexValues();
257 Long64_t *oldValues2 = GetIndexValuesMinor();
258
259 fIndex = new Long64_t[fN];
260 fIndexValues = new Long64_t[fN];
262
263 // Copy data
264 Long_t size = sizeof(Long64_t) * oldn;
265 Long_t add_size = sizeof(Long64_t) * add->GetN();
266
267 memcpy(fIndex,oldIndex, size);
268 memcpy(fIndexValues,oldValues, size);
269 memcpy(fIndexValuesMinor,oldValues2, size);
270
271 Long64_t *addIndex = ti_add->GetIndex();
272 Long64_t *addValues = ti_add->GetIndexValues();
273 Long64_t *addValues2 = ti_add->GetIndexValuesMinor();
274
275 memcpy(fIndex + oldn, addIndex, add_size);
276 memcpy(fIndexValues + oldn, addValues, add_size);
277 memcpy(fIndexValuesMinor + oldn, addValues2, add_size);
278 for(Int_t i = 0; i < add->GetN(); i++) {
279 fIndex[oldn + i] += oldn;
280 }
281
282 delete [] oldIndex;
283 delete [] oldValues;
284 delete [] oldValues2;
285 }
286
287 // Sort.
288 if (!delaySort) {
289 Long64_t *addValues = GetIndexValues();
290 Long64_t *addValues2 = GetIndexValuesMinor();
291 Long64_t *ind = fIndex;
292 Long64_t *conv = new Long64_t[fN];
293
294 for(Long64_t i = 0; i < fN; i++) { conv[i] = i; }
295 std::sort(conv, conv+fN, IndexSortComparator(addValues, addValues2) );
296 //Long64_t *w = fIndexValues;
297 //TMath::Sort(fN,w,conv,0);
298
299 fIndex = new Long64_t[fN];
300 fIndexValues = new Long64_t[fN];
302
303 for (Int_t i=0;i<fN;i++) {
304 fIndex[i] = ind[conv[i]];
305 fIndexValues[i] = addValues[conv[i]];
306 fIndexValuesMinor[i] = addValues2[conv[i]];
307 }
308 delete [] addValues;
309 delete [] addValues2;
310 delete [] ind;
311 delete [] conv;
312 }
313}
314
315
316
317////////////////////////////////////////////////////////////////////////////////
318/// conversion from old 64bit indexes
319/// return true if index was converted
320
322{
323 if( !fIndexValuesMinor && fN ) {
325 for(int i=0; i<fN; i++) {
326 fIndexValuesMinor[i] = (fIndexValues[i] & 0x7fffffff);
327 fIndexValues[i] >>= 31;
328 }
329 return true;
330 }
331 return false;
332}
333
334
335
336////////////////////////////////////////////////////////////////////////////////
337/// Returns the entry number in this (friend) Tree corresponding to entry in
338/// the master Tree 'parent'.
339/// In case this (friend) Tree and 'master' do not share an index with the same
340/// major and minor name, the entry serial number in the (friend) tree
341/// and in the master Tree are assumed to be the same
342
344{
345 if (!parent) return -3;
346 // We reached the end of the parent tree
347 Long64_t pentry = parent->GetReadEntry();
348 if (pentry >= parent->GetEntries())
349 return -2;
350 GetMajorFormulaParent(parent);
351 GetMinorFormulaParent(parent);
352 if (!fMajorFormulaParent || !fMinorFormulaParent) return -1;
354 // The Tree Index in the friend has a pair majorname,minorname
355 // not available in the parent Tree T.
356 // if the friend Tree has less entries than the parent, this is an error
357 if (pentry >= fTree->GetEntries()) return -2;
358 // otherwise we ignore the Tree Index and return the entry number
359 // in the parent Tree.
360 return pentry;
361 }
362
363 // majorname, minorname exist in the parent Tree
364 // we find the current values pair majorv,minorv in the parent Tree
367 Long64_t majorv = (Long64_t)majord;
368 Long64_t minorv = (Long64_t)minord;
369 // we check if this pair exist in the index.
370 // if yes, we return the corresponding entry number
371 // if not the function returns -1
372 return fTree->GetEntryNumberWithIndex(majorv,minorv);
373}
374
375
376////////////////////////////////////////////////////////////////////////////////
377/// find position where major|minor values are in the IndexValues tables
378/// this is the index in IndexValues table, not entry# !
379/// use lower_bound STD algorithm.
380
382{
383 Long64_t mid, step, pos = 0, count = fN;
384 // find lower bound using bisection
385 while( count > 0 ) {
386 step = count / 2;
387 mid = pos + step;
388 // check if *mid < major|minor
389 if( fIndexValues[mid] < major
390 || ( fIndexValues[mid] == major && fIndexValuesMinor[mid] < minor ) ) {
391 pos = mid+1;
392 count -= step + 1;
393 } else
394 count = step;
395 }
396 return pos;
397}
398
399
400////////////////////////////////////////////////////////////////////////////////
401/// Return entry number corresponding to major and minor number.
402/// Note that this function returns only the entry number, not the data
403/// To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
404/// the BuildIndex function has created a table of Double_t* of sorted values
405/// corresponding to val = major<<31 + minor;
406/// The function performs binary search in this sorted table.
407/// If it finds a pair that maches val, it returns directly the
408/// index in the table.
409/// If an entry corresponding to major and minor is not found, the function
410/// returns the index of the major,minor pair immediately lower than the
411/// requested value, ie it will return -1 if the pair is lower than
412/// the first entry in the index.
413///
414/// See also GetEntryNumberWithIndex
415
417{
418 if (fN == 0) return -1;
419
420 Long64_t pos = FindValues(major, minor);
421 if( pos < fN && fIndexValues[pos] == major && fIndexValuesMinor[pos] == minor )
422 return fIndex[pos];
423 if( --pos < 0 )
424 return -1;
425 return fIndex[pos];
426}
427
428
429////////////////////////////////////////////////////////////////////////////////
430/// Return entry number corresponding to major and minor number.
431/// Note that this function returns only the entry number, not the data
432/// To read the data corresponding to an entry number, use TTree::GetEntryWithIndex
433/// the BuildIndex function has created a table of Double_t* of sorted values
434/// corresponding to val = major<<31 + minor;
435/// The function performs binary search in this sorted table.
436/// If it finds a pair that maches val, it returns directly the
437/// index in the table, otherwise it returns -1.
438///
439/// See also GetEntryNumberWithBestIndex
440
442{
443 if (fN == 0) return -1;
444
445 Long64_t pos = FindValues(major, minor);
446 if( pos < fN && fIndexValues[pos] == major && fIndexValuesMinor[pos] == minor )
447 return fIndex[pos];
448 return -1;
449}
450
451
452////////////////////////////////////////////////////////////////////////////////
453
455{
456 return fIndexValuesMinor;
457}
458
459
460
461////////////////////////////////////////////////////////////////////////////////
462/// Return a pointer to the TreeFormula corresponding to the majorname.
463
465{
466 if (!fMajorFormula) {
469 }
470 return fMajorFormula;
471}
472
473////////////////////////////////////////////////////////////////////////////////
474/// Return a pointer to the TreeFormula corresponding to the minorname.
475
477{
478 if (!fMinorFormula) {
481 }
482 return fMinorFormula;
483}
484
485////////////////////////////////////////////////////////////////////////////////
486/// Return a pointer to the TreeFormula corresponding to the majorname in parent tree.
487
489{
490 if (!fMajorFormulaParent) {
491 // Prevent TTreeFormula from finding any of the branches in our TTree even if it
492 // is a friend of the parent TTree.
494 fMajorFormulaParent = new TTreeFormula("MajorP",fMajorName.Data(),const_cast<TTree*>(parent));
496 }
497 if (fMajorFormulaParent->GetTree() != parent) {
498 fMajorFormulaParent->SetTree(const_cast<TTree*>(parent));
500 }
501 return fMajorFormulaParent;
502}
503
504////////////////////////////////////////////////////////////////////////////////
505/// Return a pointer to the TreeFormula corresponding to the minorname in parent tree.
506
508{
509 if (!fMinorFormulaParent) {
510 // Prevent TTreeFormula from finding any of the branches in our TTree even if it
511 // is a friend of the parent TTree.
513 fMinorFormulaParent = new TTreeFormula("MinorP",fMinorName.Data(),const_cast<TTree*>(parent));
515 }
516 if (fMinorFormulaParent->GetTree() != parent) {
517 fMinorFormulaParent->SetTree(const_cast<TTree*>(parent));
519 }
520 return fMinorFormulaParent;
521}
522
523////////////////////////////////////////////////////////////////////////////////
524/// Return true if index can be applied to the TTree
525
526bool TTreeIndex::IsValidFor(const TTree *parent)
527{
528 auto *majorFormula = GetMajorFormulaParent(parent);
529 auto *minorFormula = GetMinorFormulaParent(parent);
530 if ((majorFormula == nullptr || majorFormula->GetNdim() == 0) ||
531 (minorFormula == nullptr || minorFormula->GetNdim() == 0))
532 return false;
533 return true;
534}
535
536////////////////////////////////////////////////////////////////////////////////
537/// Print the table with : serial number, majorname, minorname.
538/// - if option = "10" print only the first 10 entries
539/// - if option = "100" print only the first 100 entries
540/// - if option = "1000" print only the first 1000 entries
541
543{
544 TString opt = option;
545 bool printEntry = false;
546 Long64_t n = fN;
547 if (opt.Contains("10")) n = 10;
548 if (opt.Contains("100")) n = 100;
549 if (opt.Contains("1000")) n = 1000;
550 if (opt.Contains("all")) {
551 printEntry = true;
552 }
553
554 if (printEntry) {
555 Printf("\n*****************************************************************");
556 Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
557 Printf("*****************************************************************");
558 Printf("%8s : %16s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data(),"entry number");
559 Printf("*****************************************************************");
560 for (Long64_t i=0;i<n;i++) {
561 Printf("%8lld : %8lld : %8lld : %8lld",
562 i, fIndexValues[i], GetIndexValuesMinor()[i], fIndex[i]);
563 }
564
565 } else {
566 Printf("\n**********************************************");
567 Printf("* Index of Tree: %s/%s",fTree->GetName(),fTree->GetTitle());
568 Printf("**********************************************");
569 Printf("%8s : %16s : %16s","serial",fMajorName.Data(),fMinorName.Data());
570 Printf("**********************************************");
571 for (Long64_t i=0;i<n;i++) {
572 Printf("%8lld : %8lld : %8lld",
574 }
575 }
576}
577
578////////////////////////////////////////////////////////////////////////////////
579/// Stream an object of class TTreeIndex.
580/// Note that this Streamer should be changed to an automatic Streamer
581/// once TStreamerInfo supports an index of type Long64_t
582
584{
585 UInt_t R__s, R__c;
586 if (R__b.IsReading()) {
587 Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
589 fMajorName.Streamer(R__b);
590 fMinorName.Streamer(R__b);
591 R__b >> fN;
592 fIndexValues = new Long64_t[fN];
594 if( R__v > 1 ) {
597 } else {
599 }
600 fIndex = new Long64_t[fN];
601 R__b.ReadFastArray(fIndex,fN);
602 R__b.CheckByteCount(R__s, R__c, TTreeIndex::IsA());
603 } else {
604 R__c = R__b.WriteVersion(TTreeIndex::IsA(), true);
606 fMajorName.Streamer(R__b);
607 fMinorName.Streamer(R__b);
608 R__b << fN;
611 R__b.WriteFastArray(fIndex, fN);
612 R__b.SetByteCount(R__c, true);
613 }
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Called by TChain::LoadTree when the parent chain changes it's tree.
618
620{
624 if (parent) fMajorFormulaParent->SetTree(const_cast<TTree*>(parent));
626 }
628 if (parent) fMinorFormulaParent->SetTree(const_cast<TTree*>(parent));
630 }
631}
632////////////////////////////////////////////////////////////////////////////////
633/// this function is called by TChain::LoadTree and TTreePlayer::UpdateFormulaLeaves
634/// when a new Tree is loaded.
635/// Because Trees in a TChain may have a different list of leaves, one
636/// must update the leaves numbers in the TTreeFormula used by the TreeIndex.
637
639{
640 fTree = T;
641}
642
643////////////////////////////////////////////////////////////////////////////////
644/// \brief Create a deep copy of the TTreeIndex
645/// \param[in] newname A new name for the index
646///
647/// The new index is allocated on the heap without being managed. Also, it is
648/// not attached to any tree. It is the responsibility of the caller to manage
649/// its lifetime and attach it to a tree if necessary.
650TObject *TTreeIndex::Clone(const char *newname) const
651{
652 auto index = new TTreeIndex();
653 index->SetName(newname && std::strlen(newname) ? newname : GetName());
654 index->SetTitle(GetTitle());
655
656 // Note that the TTreeFormula * data members are not cloned since they would
657 // need the attached tree data member to function properly.
658 index->fMajorName = fMajorName;
659 index->fMinorName = fMinorName;
660
661 if (fN == 0)
662 return index;
663
664 index->fN = fN;
665
666 index->fIndexValues = new Long64_t[index->fN];
667 std::copy(fIndexValues, fIndexValues + fN, index->fIndexValues);
668
669 index->fIndexValuesMinor = new Long64_t[index->fN];
670 std::copy(fIndexValuesMinor, fIndexValuesMinor + fN, index->fIndexValuesMinor);
671
672 index->fIndex = new Long64_t[index->fN];
673 std::copy(fIndex, fIndex + fN, index->fIndex);
674
675 return index;
676}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
short Version_t
Definition RtypesCore.h:65
long Long_t
Definition RtypesCore.h:54
long double LongDouble_t
Definition RtypesCore.h:61
long long Long64_t
Definition RtypesCore.h:80
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
void Printf(const char *fmt,...)
Formats a string in a circular formatting buffer and prints the string.
Definition TString.cxx:2503
virtual Int_t GetNdim() const
Definition TFormula.h:237
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Version_t ReadVersion(UInt_t *start=nullptr, UInt_t *bcnt=nullptr, const TClass *cl=nullptr)=0
virtual void SetByteCount(UInt_t cntpos, Bool_t packInVersion=kFALSE)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
virtual void ReadFastArray(Bool_t *b, Int_t n)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual void WriteFastArray(const Bool_t *b, Long64_t n)=0
virtual UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE)=0
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:973
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
void MakeZombie()
Definition TObject.h:53
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
virtual void Streamer(TBuffer &)
Stream a string object.
Definition TString.cxx:1412
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
Used to pass a selection expression to the Tree drawing routine.
virtual void SetTree(TTree *tree)
T EvalInstance(Int_t i=0, const char *stringStack[]=nullptr)
Evaluate this treeformula.
void SetQuickLoad(bool quick)
virtual void UpdateFormulaLeaves()
This function is called TTreePlayer::UpdateFormulaLeaves, itself called by TChain::LoadTree when a ne...
virtual TTree * GetTree() const
A Tree Index with majorname and minorname.
Definition TTreeIndex.h:29
TTreeIndex()
Default constructor for TTreeIndex.
virtual Long64_t * GetIndexValues() const
Definition TTreeIndex.h:60
virtual Long64_t * GetIndexValuesMinor() const
TTreeFormula * fMajorFormula
! Pointer to major TreeFormula
Definition TTreeIndex.h:37
TTreeFormula * fMajorFormulaParent
! Pointer to major TreeFormula in Parent tree (if any)
Definition TTreeIndex.h:39
Long64_t * fIndex
[fN] Index of sorted values
Definition TTreeIndex.h:36
Long64_t GetEntryNumberWithIndex(Long64_t major, Long64_t minor) const override
Return entry number corresponding to major and minor number.
TTreeFormula * GetMajorFormulaParent(const TTree *parent)
Return a pointer to the TreeFormula corresponding to the majorname in parent tree.
void SetTree(TTree *T) override
this function is called by TChain::LoadTree and TTreePlayer::UpdateFormulaLeaves when a new Tree is l...
bool IsValidFor(const TTree *parent) override
Return true if index can be applied to the TTree.
Long64_t fN
Number of entries.
Definition TTreeIndex.h:33
TClass * IsA() const override
Definition TTreeIndex.h:73
virtual Long64_t * GetIndex() const
Definition TTreeIndex.h:59
bool ConvertOldToNew()
conversion from old 64bit indexes return true if index was converted
TTreeFormula * fMinorFormula
! Pointer to minor TreeFormula
Definition TTreeIndex.h:38
void Append(const TVirtualIndex *, bool delaySort=false) override
Append 'add' to this index.
void UpdateFormulaLeaves(const TTree *parent) override
Called by TChain::LoadTree when the parent chain changes it's tree.
virtual TTreeFormula * GetMajorFormula()
Return a pointer to the TreeFormula corresponding to the majorname.
TTreeFormula * GetMinorFormulaParent(const TTree *parent)
Return a pointer to the TreeFormula corresponding to the minorname in parent tree.
Long64_t GetEntryNumberWithBestIndex(Long64_t major, Long64_t minor) const override
Return entry number corresponding to major and minor number.
TObject * Clone(const char *newname="") const override
Create a deep copy of the TTreeIndex.
Long64_t GetEntryNumberFriend(const TTree *parent) override
Returns the entry number in this (friend) Tree corresponding to entry in the master Tree 'parent'.
TString fMinorName
Index minor name.
Definition TTreeIndex.h:32
void Print(Option_t *option="") const override
Print the table with : serial number, majorname, minorname.
void Streamer(TBuffer &) override
Stream an object of class TTreeIndex.
virtual TTreeFormula * GetMinorFormula()
Return a pointer to the TreeFormula corresponding to the minorname.
Long64_t * fIndexValues
[fN] Sorted index values, higher 64bits
Definition TTreeIndex.h:34
TString fMajorName
Index major name.
Definition TTreeIndex.h:31
Long64_t * fIndexValuesMinor
[fN] Sorted index values, lower 64bits
Definition TTreeIndex.h:35
~TTreeIndex() override
Destructor.
TTreeFormula * fMinorFormulaParent
! Pointer to minor TreeFormula in Parent tree (if any)
Definition TTreeIndex.h:40
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...
Helper class to prevent infinite recursion in the usage of TTree Friends.
Definition TTree.h:188
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Long64_t GetEntryNumberWithIndex(Long64_t major, Long64_t minor=0) const
Return entry number corresponding to major and minor number.
Definition TTree.cxx:5910
virtual TVirtualIndex * GetTreeIndex() const
Definition TTree.h:518
virtual Long64_t GetEntries() const
Definition TTree.h:463
virtual Long64_t GetReadEntry() const
Definition TTree.h:509
virtual Long64_t LoadTree(Long64_t entry)
Set current entry.
Definition TTree.cxx:6473
virtual Int_t GetTreeNumber() const
Definition TTree.h:519
@ kFindBranch
Definition TTree.h:212
@ kFindLeaf
Definition TTree.h:213
@ kGetBranch
Definition TTree.h:215
@ kGetLeaf
Definition TTree.h:220
virtual void SetTreeIndex(TVirtualIndex *index)
The current TreeIndex is replaced by the new index.
Definition TTree.cxx:9346
Abstract interface for Tree Index.
void Streamer(TBuffer &) override
Stream an object of class TObject.
virtual Long64_t GetN() const =0
TClass * IsA() const override
const Int_t n
Definition legend1.C:16
bool operator()(Index i1, Index i2)
IndexSortComparator(Long64_t *major, Long64_t *minor)