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