Logo ROOT   6.07/09
Reference Guide
TFormLeafInfo.cxx
Go to the documentation of this file.
1 // @(#)root/treeplayer:$Id$
2 // Author: Philippe Canal 01/06/2004
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers and al. *
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 TFormLeafInfo
13 This class is a small helper class to implement reading a data member
14 on an object stored in a TTree.
15 
16 TTreeFormula now relies on a variety of TFormLeafInfo classes to handle the
17 reading of the information. Here is the list of theses classes:
18  - TFormLeafInfo
19  - TFormLeafInfoDirect
20  - TFormLeafInfoNumerical
21  - TFormLeafInfoClones
22  - TFormLeafInfoCollection
23  - TFormLeafInfoPointer
24  - TFormLeafInfoMethod
25  - TFormLeafInfoMultiVarDim
26  - TFormLeafInfoMultiVarDimDirect
27  - TFormLeafInfoCast
28 
29 The following method are available from the TFormLeafInfo interface:
30 
31  - AddOffset(Int_t offset, TStreamerElement* element)
32  - GetCounterValue(TLeaf* leaf) : return the size of the array pointed to.
33  - GetObjectAddress(TLeafElement* leaf) : Returns the the location of the object pointed to.
34  - GetMultiplicity() : Returns info on the variability of the number of elements
35  - GetNdata(TLeaf* leaf) : Returns the number of elements
36  - GetNdata() : Used by GetNdata(TLeaf* leaf)
37  - GetValue(TLeaf *leaf, Int_t instance = 0) : Return the value
38  - GetValuePointer(TLeaf *leaf, Int_t instance = 0) : Returns the address of the value
39  - GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0) : Returns the address of the value of 'this' LeafInfo
40  - IsString()
41  - ReadValue(char *where, Int_t instance = 0) : Internal function to interpret the location 'where'
42  - Update() : react to the possible loading of a shared library.
43 */
44 
45 #include "TFormLeafInfo.h"
46 
47 #include "TROOT.h"
48 #include "TArrayI.h"
49 #include "TClonesArray.h"
50 #include "TError.h"
51 #include "TInterpreter.h"
52 #include "TLeafObject.h"
53 #include "TMethod.h"
54 #include "TMethodCall.h"
55 #include "TTree.h"
57 
58 
59 #define INSTANTIATE_READVAL(CLASS) \
60  template Double_t CLASS::ReadValueImpl<Double_t>(char*, Int_t); \
61  template Long64_t CLASS::ReadValueImpl<Long64_t>(char*, Int_t); \
62  template LongDouble_t CLASS::ReadValueImpl<LongDouble_t>(char*, Int_t) // no semicolon
63 
64 
65 #define INSTANTIATE_GETVAL(CLASS) \
66  template Double_t CLASS::GetValueImpl<Double_t>(TLeaf*, Int_t); \
67  template Long64_t CLASS::GetValueImpl<Long64_t>(TLeaf*, Int_t); \
68  template LongDouble_t CLASS::GetValueImpl<LongDouble_t>(TLeaf*, Int_t) // no semicolon
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 /// Constructor.
72 
74  TStreamerElement* element) :
75  fClass(classptr),fOffset(offset),fElement(element),
76  fCounter(0), fNext(0),fMultiplicity(0)
77 {
78  if (fClass) fClassName = fClass->GetName();
79  if (fElement) {
81  }
82 }
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 ///Constructor.
86 
88 {
89  // Deep copy the pointers.
90  if (orig.fCounter) fCounter = orig.fCounter->DeepCopy();
91  if (orig.fNext) fNext = orig.fNext->DeepCopy();
92 }
93 
94 ////////////////////////////////////////////////////////////////////////////////
95 /// Exception safe assignment operator.
96 
98 {
99  TFormLeafInfo tmp(other);
100  Swap(tmp);
101  return *this;
102 }
103 
104 ////////////////////////////////////////////////////////////////////////////////
105 
107 {
108  std::swap(fClass,other.fClass);
109  std::swap(fOffset,other.fOffset);
110  std::swap(fElement,other.fElement);
111  std::swap(fCounter,other.fCounter);
112  std::swap(fNext,other.fNext);
113  TString tmp(fClassName);
114  fClassName = other.fClassName;
115  other.fClassName = tmp;
116 
117  tmp = fElementName;
118  fElementName = other.fElementName;
119  other.fElementName = tmp;
121 }
122 
123 ////////////////////////////////////////////////////////////////////////////////
124 /// Make a complete copy of this FormLeafInfo and all its content.
125 
127 {
128  return new TFormLeafInfo(*this);
129 }
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// Delete this object and all its content
133 
135 {
136  delete fCounter;
137  delete fNext;
138 }
139 
140 ////////////////////////////////////////////////////////////////////////////////
141 /// Increase the offset of this element. This intended to be the offset
142 /// from the start of the object to which the data member belongs.
143 
145 {
146  fOffset += offset;
147  fElement = element;
148  if (fElement ) {
149  // fElementClassOwnerName = cl->GetName();
150  fElementName.Append(".").Append(element->GetName());
151  }
152 }
153 
154 ////////////////////////////////////////////////////////////////////////////////
155 /// Return the current length of the array.
156 
158 {
159  Int_t len = 1;
160  if (fNext) len = fNext->GetArrayLength();
161  if (fElement) {
162  Int_t elen = fElement->GetArrayLength();
163  if (elen || fElement->IsA() == TStreamerBasicPointer::Class() )
164  len *= fElement->GetArrayLength();
165  }
166  return len;
167 }
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 /// Get the class of the underlying data.
171 
173 {
174  if (fNext) return fNext->GetClass();
175  if (fElement) return fElement->GetClassPointer();
176  return fClass;
177 }
178 
179 ////////////////////////////////////////////////////////////////////////////////
180 /// Returns the the location of the object pointed to.
181 /// Modify instance if the object is part of an array.
182 
184 {
185  TBranchElement* branch = (TBranchElement*) leaf->GetBranch();
186  Int_t id = branch->GetID();
187  if (id < 0) {
188  // Branch is a top-level branch.
189  if (branch->GetTree()->GetMakeClass()) {
190  // Branch belongs to a MakeClass tree.
191  return branch->GetAddress();
192  } else {
193  return branch->GetObject();
194  }
195  }
196  TStreamerInfo* info = branch->GetInfo();
197  Int_t offset = 0;
198  if (id > -1) {
199  // Branch is *not* a top-level branch.
200  offset = info->TStreamerInfo::GetElementOffset(id);
201  }
202  char* address = 0;
203  // Branch is *not* a top-level branch.
204  if (branch->GetTree()->GetMakeClass()) {
205  // Branch belongs to a MakeClass tree.
206  address = (char*) branch->GetAddress();
207  } else {
208  address = (char*) branch->GetObject();
209  }
210  char* thisobj = 0;
211  if (!address) {
212  // FIXME: This makes no sense, if the branch address is not set, then object will not be set either.
213  thisobj = branch->GetObject();
214  } else {
215  Int_t type = -1;
216  if (id > -1) {
217  // Note this is somewhat slow
218  type = info->TStreamerInfo::GetElement(id)->GetNewType();
219  }
220  switch (type) {
226  Error("GetValuePointer", "Type (%d) not yet supported\n", type);
227  break;
228 
232  {
233  // An array of objects.
234  Int_t index;
235  Int_t sub_instance;
236  Int_t len = GetArrayLength();
237  if (len) {
238  index = instance / len;
239  sub_instance = instance % len;
240  } else {
241  index = instance;
242  sub_instance = 0;
243  }
244  thisobj = address + offset + (index * fClass->Size());
245  instance = sub_instance;
246  break;
247  }
248 
254  case TStreamerInfo::kAny:
255  case TStreamerInfo::kSTL:
256  // A single object.
257  thisobj = address + offset;
258  break;
259 
265  case TStreamerInfo::kInt:
293  // A simple type, or an array of a simple type.
294  thisobj = address + offset;
295  break;
296 
297  default:
298  // Everything else is a pointer to something.
299  thisobj = *((char**) (address + offset));
300  break;
301  }
302  }
303  return thisobj;
304 }
305 
306 ////////////////////////////////////////////////////////////////////////////////
307 /// Reminder of the meaning of fMultiplicity:
308 /// - -1: Only one or 0 element per entry but contains variable length array!
309 /// - 0: Only one element per entry, no variable length array
310 /// - 1: loop over the elements of a variable length array
311 /// - 2: loop over elements of fixed length array (nData is the same for all entry)
312 
314 {
315  // Currently only TFormLeafInfoCast uses this field.
316  return fMultiplicity;
317 }
318 
320 {
321  // Get the number of element in the entry.
322 
323  GetCounterValue(leaf);
324  GetValue(leaf);
325  return GetNdata();
326 }
327 
328 ////////////////////////////////////////////////////////////////////////////////
329 /// Get the number of element in the entry.
330 
332 {
333  if (fNext) return fNext->GetNdata();
334  return 1;
335 }
336 
337 ////////////////////////////////////////////////////////////////////////////////
338 /// Return true if any of underlying data has a array size counter
339 
341 {
342  Bool_t result = kFALSE;
343  if (fNext) result = fNext->HasCounter();
344  return fCounter!=0 || result;
345 }
346 
347 ////////////////////////////////////////////////////////////////////////////////
348 /// Return true if the underlying data is a string
349 
351 {
352  if (fNext) return fNext->IsString();
353  if (!fElement) return kFALSE;
354 
355  switch (fElement->GetNewType()) {
356  // basic types
357  case kChar_t:
358  // This is new in ROOT 3.02/05
359  return kFALSE;
361  // This is new in ROOT 3.02/05
362  return kTRUE;
364  return kTRUE;
365  default:
366  return kFALSE;
367  }
368 }
369 
370 ////////////////////////////////////////////////////////////////////////////////
371 /// Return true if the underlying data is an integral value
372 
374 {
375  if (fNext) return fNext->IsInteger();
376  if (!fElement) return kFALSE;
377 
378  Int_t atype = fElement->GetNewType();
379  if (TStreamerInfo::kOffsetL < atype &&
380  atype < TStreamerInfo::kOffsetP ) {
381  atype -= TStreamerInfo::kOffsetL;
382  } else if (TStreamerInfo::kOffsetP < atype &&
383  atype < TStreamerInfo::kObject) {
384  atype -= TStreamerInfo::kOffsetP;
385  }
386 
387  switch (atype) {
388  // basic types
395  case TStreamerInfo::kInt:
401  return kTRUE;
403  return kTRUE; // For consistency with the leaf list method and proper axis setting
408  return kFALSE;
409  default:
410  return kFALSE;
411  }
412 }
413 
414 ////////////////////////////////////////////////////////////////////////////////
415 /// Method for multiple variable dimensions.
416 
418 {
419  if (fNext) return fNext->GetPrimaryIndex();
420  return -1;
421 }
422 
423 ////////////////////////////////////////////////////////////////////////////////
424 /// Return the index of the dimension which varies
425 /// for each elements of an enclosing array (typically a TClonesArray)
426 
428 {
429  if (fNext) return fNext->GetVarDim();
430  else return -1;
431 }
432 
433 ////////////////////////////////////////////////////////////////////////////////
434 /// Return the virtual index (for this expression) of the dimension which varies
435 /// for each elements of an enclosing array (typically a TClonesArray)
436 
438 {
439  if (fNext) return fNext->GetVirtVarDim();
440  else return -1;
441 }
442 
443 ////////////////////////////////////////////////////////////////////////////////
444 /// For the current entry, and the value 'index' for the main array,
445 /// return the size of the secondary variable dimension of the 'array'.
446 
448 {
449  if (fNext) return fNext->GetSize(index);
450  else return 0;
451 }
452 
453 ////////////////////////////////////////////////////////////////////////////////
454 /// Total all the elements that are available for the current entry
455 /// for the secondary variable dimension.
456 
458 {
459  if (fNext) return fNext->GetSumOfSizes();
460  else return 0;
461 }
462 
463 ////////////////////////////////////////////////////////////////////////////////
464 /// Load the current array sizes
465 
467 {
468  if (fNext) fNext->LoadSizes(branch);
469 }
470 
471 ////////////////////////////////////////////////////////////////////////////////
472 /// Set the primary index value
473 
475 {
476  if (fNext) fNext->SetPrimaryIndex(index);
477 }
478 
479 ////////////////////////////////////////////////////////////////////////////////
480 /// Set the primary index value
481 
483 {
484  if (fNext) fNext->SetSecondaryIndex(index);
485 }
486 
487 ////////////////////////////////////////////////////////////////////////////////
488 /// Set the current size of the arrays
489 
491 {
492  if (fNext) fNext->SetSize(index, val);
493 }
494 
495 ////////////////////////////////////////////////////////////////////////////////
496 /// Set the current sizes of the arrays
497 
499 {
500  if (fNext) fNext->UpdateSizes(garr);
501 }
502 
503 ////////////////////////////////////////////////////////////////////////////////
504 /// We reloading all cached information in case the underlying class
505 /// information has changed (for example when changing from the 'emulated'
506 /// class to the real class.
507 
509 {
510  if (fClass) {
511  TClass * new_class = TClass::GetClass(fClassName);
512  if (new_class==fClass) {
513  if (fNext) fNext->Update();
514  if (fCounter) fCounter->Update();
515  return kFALSE;
516  }
517  fClass = new_class;
518  }
519  if (fElement && fClass) {
520  TClass *cl = fClass;
521  // We have to drill down the element name within the class.
522  Int_t offset,i;
523  TStreamerElement* element;
524  char * current;
525  Int_t nchname = fElementName.Length();
526  char * work = new char[nchname+2];
527  for (i=0, current = &(work[0]), fOffset=0; i<nchname+1;i++ ) {
528  if (i==nchname || fElementName[i]=='.') {
529  // A delimiter happened let's see if what we have seen
530  // so far does point to a data member.
531  *current = '\0';
532  element = ((TStreamerInfo*)cl->GetStreamerInfo())->GetStreamerElement(work,offset);
533  if (element) {
534  Int_t type = element->GetNewType();
535  if (type<60) {
536  fOffset += offset;
537  } else if (type == TStreamerInfo::kBase ||
538  type == TStreamerInfo::kAny ||
539  type == TStreamerInfo::kObject ||
540  type == TStreamerInfo::kTString ||
541  type == TStreamerInfo::kTNamed ||
542  type == TStreamerInfo::kTObject ||
543  type == TStreamerInfo::kObjectp ||
544  type == TStreamerInfo::kObjectP ||
547  type == TStreamerInfo::kAnyp ||
548  type == TStreamerInfo::kAnyP ||
552  type == TStreamerInfo::kSTL ||
553  type == TStreamerInfo::kSTLp ) {
554  fOffset += offset;
555  cl = element->GetClassPointer();
556  }
557  fElement = element;
558  current = &(work[0]);
559  }
560  } else {
561  if (i<nchname) *current++ = fElementName[i];
562  }
563  }
564  delete [] work;
565  }
566  if (fNext) fNext->Update();
567  if (fCounter) fCounter->Update();
568  return kTRUE;
569 }
570 
571 ////////////////////////////////////////////////////////////////////////////////
572 /// Return the size of the underlying array for the current entry in the TTree.
573 
575  if (!fCounter) {
576  if (fNext && fNext->HasCounter()) {
577  char *where = (char*)GetLocalValuePointer(leaf,0);
578  return fNext->ReadCounterValue(where);
579  } else return 1;
580  }
581  return (Int_t)fCounter->GetValue(leaf);
582 }
583 
584 ////////////////////////////////////////////////////////////////////////////////
585 /// Return the size of the underlying array for the current entry in the TTree.
586 
588 {
589  if (!fCounter) {
590  if (fNext) {
591  char *next = (char*)GetLocalValuePointer(where,0);
592  return fNext->ReadCounterValue(next);
593  } else return 1;
594  }
595  return (Int_t)fCounter->ReadValue(where,0);
596 }
597 
598 ////////////////////////////////////////////////////////////////////////////////
599 /// returns the address of the value pointed to by the
600 /// TFormLeafInfo.
601 
603 {
604  char *thisobj = 0;
605  if (leaf->InheritsFrom(TLeafObject::Class()) ) {
606  thisobj = (char*)((TLeafObject*)leaf)->GetObject();
607  } else {
608  thisobj = GetObjectAddress((TLeafElement*)leaf, instance); // instance might be modified
609  }
610  if (!thisobj) return 0;
611  return GetLocalValuePointer(thisobj, instance);
612 }
613 
614 ////////////////////////////////////////////////////////////////////////////////
615 /// returns the address of the value pointed to by the
616 /// serie of TFormLeafInfo.
617 
619 {
620  char *thisobj = (char*)GetLocalValuePointer(leaf,instance);
621  if (fNext) return fNext->GetValuePointer(thisobj,instance);
622  else return thisobj;
623 }
624 
625 ////////////////////////////////////////////////////////////////////////////////
626 /// returns the address of the value pointed to by the
627 /// TFormLeafInfo.
628 
629 void* TFormLeafInfo::GetValuePointer(char *thisobj, Int_t instance)
630 {
631  char *where = (char*)GetLocalValuePointer(thisobj,instance);
632  if (fNext) return fNext->GetValuePointer(where,instance);
633  else return where;
634 }
635 
636 ////////////////////////////////////////////////////////////////////////////////
637 /// returns the address of the value pointed to by the
638 /// TFormLeafInfo.
639 
640 void* TFormLeafInfo::GetLocalValuePointer(char *thisobj, Int_t instance)
641 {
642  if (fElement==0 || thisobj==0) return thisobj;
643 
644  switch (fElement->GetNewType()) {
645  // basic types
651  case TStreamerInfo::kInt:
663  return (Int_t*)(thisobj+fOffset);
664 
665  // array of basic types array[8]
667  {Bool_t *val = (Bool_t*)(thisobj+fOffset); return &(val[instance]);}
669  {Char_t *val = (Char_t*)(thisobj+fOffset); return &(val[instance]);}
671  {Short_t *val = (Short_t*)(thisobj+fOffset); return &(val[instance]);}
673  {Int_t *val = (Int_t*)(thisobj+fOffset); return &(val[instance]);}
675  {Long_t *val = (Long_t*)(thisobj+fOffset); return &(val[instance]);}
677  {Long64_t *val = (Long64_t*)(thisobj+fOffset); return &(val[instance]);}
679  {Float_t *val = (Float_t*)(thisobj+fOffset); return &(val[instance]);}
681  {Float_t *val = (Float_t*)(thisobj+fOffset); return &(val[instance]);}
683  {Double_t *val = (Double_t*)(thisobj+fOffset); return &(val[instance]);}
685  {Double_t *val = (Double_t*)(thisobj+fOffset); return &(val[instance]);}
687  {UChar_t *val = (UChar_t*)(thisobj+fOffset); return &(val[instance]);}
689  {UShort_t *val = (UShort_t*)(thisobj+fOffset); return &(val[instance]);}
691  {UInt_t *val = (UInt_t*)(thisobj+fOffset); return &(val[instance]);}
693  {ULong_t *val = (ULong_t*)(thisobj+fOffset); return &(val[instance]);}
695  {ULong64_t *val = (ULong64_t*)(thisobj+fOffset); return &(val[instance]);}
696 
697 #define GET_ARRAY(TYPE_t) \
698  { \
699  Int_t len, sub_instance, index; \
700  if (fNext) len = fNext->GetArrayLength(); \
701  else len = 1; \
702  if (len) { \
703  index = instance / len; \
704  sub_instance = instance % len; \
705  } else { \
706  index = instance; \
707  sub_instance = 0; \
708  } \
709  TYPE_t **val = (TYPE_t**)(thisobj+fOffset); \
710  return &((val[sub_instance])[index]); \
711  }
712 
713  // pointer to an array of basic types array[n]
729 
731  {char **stringp = (char**)(thisobj+fOffset); return *stringp;}
732 
738  {TObject **obj = (TObject**)(thisobj+fOffset); return *obj; }
739 
744  case TStreamerInfo::kAny:
746  case TStreamerInfo::kSTL:
747  {TObject *obj = (TObject*)(thisobj+fOffset); return obj; }
748 
752  char *loc = thisobj+fOffset;
753 
754  Int_t len, index;
755  //Int_t sub_instance;
756 
757  if (fNext) len = fNext->GetArrayLength();
758  else len = 1;
759  if (len) {
760  index = instance / len;
761  // sub_instance = instance % len;
762  } else {
763  index = instance;
764  // sub_instance = 0;
765  }
766 
767  loc += index*fElement->GetClassPointer()->Size();
768 
769  TObject *obj = (TObject*)(loc);
770  return obj;
771  }
772 
778  {TObject *obj = (TObject*)(thisobj+fOffset); return obj; }
779 
780  case kOther_t:
781  default: return 0;
782  }
783 }
784 
785 
786 ////////////////////////////////////////////////////////////////////////////////
787 /// Return result of a leafobject method.
788 
789 template <typename T>
790 T TFormLeafInfo::GetValueImpl(TLeaf *leaf, Int_t instance)
791 {
792  char *thisobj = 0;
793  if (leaf->InheritsFrom(TLeafObject::Class()) ) {
794  thisobj = (char*)((TLeafObject*)leaf)->GetObject();
795  } else {
796  thisobj = GetObjectAddress((TLeafElement*)leaf, instance); // instance might be modified
797  }
798  if (thisobj==0) return 0;
799  return ReadTypedValue<T>(thisobj,instance);
800 }
801 
804 
805 ////////////////////////////////////////////////////////////////////////////////
806 /// Read the value at the given memory location
807 
808 template <typename T>
809 T TFormLeafInfo::ReadValueImpl(char *thisobj, Int_t instance)
810 {
811  if ( !thisobj ) {
812  Error("ReadValue","Invalid data address: result will be wrong");
813  return 0.0; // Or should throw exception/print error ?
814  }
815  if (fNext) {
816  char *nextobj = thisobj+fOffset;
817  Int_t sub_instance = instance;
822  Int_t index;
823  Int_t len = fNext->GetArrayLength();
824  if (len) {
825  index = instance / len;
826  sub_instance = instance % len;
827  } else {
828  index = instance;
829  sub_instance = 0;
830  }
831  nextobj += index*fElement->GetClassPointer()->Size();
832  }
833  return fNext->ReadTypedValue<T>(nextobj,sub_instance);
834  }
835  // return fInfo->ReadValue(thisobj+fOffset,fElement->GetNewType(),instance,1);
836  switch (fElement->GetNewType()) {
837  // basic types
838  case TStreamerInfo::kBool: return (T)(*(Bool_t*)(thisobj+fOffset));
839  case TStreamerInfo::kChar: return (T)(*(Char_t*)(thisobj+fOffset));
840  case TStreamerInfo::kUChar: return (T)(*(UChar_t*)(thisobj+fOffset));
841  case TStreamerInfo::kShort: return (T)(*(Short_t*)(thisobj+fOffset));
842  case TStreamerInfo::kUShort: return (T)(*(UShort_t*)(thisobj+fOffset));
843  case TStreamerInfo::kInt: return (T)(*(Int_t*)(thisobj+fOffset));
844  case TStreamerInfo::kUInt: return (T)(*(UInt_t*)(thisobj+fOffset));
845  case TStreamerInfo::kLong: return (T)(*(Long_t*)(thisobj+fOffset));
846  case TStreamerInfo::kULong: return (T)(*(ULong_t*)(thisobj+fOffset));
847  case TStreamerInfo::kLong64: return (T)(*(Long64_t*)(thisobj+fOffset));
848  case TStreamerInfo::kULong64: return (T)(*(Long64_t*)(thisobj+fOffset)); //cannot cast to ULong64_t with VC++6
849  case TStreamerInfo::kFloat: return (T)(*(Float_t*)(thisobj+fOffset));
850  case TStreamerInfo::kFloat16: return (T)(*(Float_t*)(thisobj+fOffset));
851  case TStreamerInfo::kDouble: return (T)(*(Double_t*)(thisobj+fOffset));
852  case TStreamerInfo::kDouble32: return (T)(*(Double_t*)(thisobj+fOffset));
853  case TStreamerInfo::kLegacyChar: return (T)(*(char*)(thisobj+fOffset));
854  case TStreamerInfo::kCounter: return (T)(*(Int_t*)(thisobj+fOffset));
855 
856  // array of basic types array[8]
858  {Bool_t *val = (Bool_t*)(thisobj+fOffset); return T(val[instance]);}
860  {Char_t *val = (Char_t*)(thisobj+fOffset); return T(val[instance]);}
862  {Short_t *val = (Short_t*)(thisobj+fOffset); return T(val[instance]);}
864  {Int_t *val = (Int_t*)(thisobj+fOffset); return T(val[instance]);}
866  {Long_t *val = (Long_t*)(thisobj+fOffset); return T(val[instance]);}
868  {Long64_t *val = (Long64_t*)(thisobj+fOffset); return T(val[instance]);}
870  {Float_t *val = (Float_t*)(thisobj+fOffset); return T(val[instance]);}
872  {Float_t *val = (Float_t*)(thisobj+fOffset); return T(val[instance]);}
874  {Double_t *val = (Double_t*)(thisobj+fOffset); return T(val[instance]);}
876  {Double_t *val = (Double_t*)(thisobj+fOffset); return T(val[instance]);}
878  {UChar_t *val = (UChar_t*)(thisobj+fOffset); return T(val[instance]);}
880  {UShort_t *val = (UShort_t*)(thisobj+fOffset); return T(val[instance]);}
882  {UInt_t *val = (UInt_t*)(thisobj+fOffset); return T(val[instance]);}
884  {ULong_t *val = (ULong_t*)(thisobj+fOffset); return T(val[instance]);}
885 #if defined(_MSC_VER) && (_MSC_VER <= 1200)
887  {Long64_t *val = (Long64_t*)(thisobj+fOffset); return T(val[instance]);}
888 #else
890  {ULong64_t *val = (ULong64_t*)(thisobj+fOffset); return T(val[instance]);}
891 #endif
892 
893 #define READ_ARRAY(TYPE_t) \
894  { \
895  Int_t len, sub_instance, index; \
896  len = GetArrayLength(); \
897  if (len) { \
898  index = instance / len; \
899  sub_instance = instance % len; \
900  } else { \
901  index = instance; \
902  sub_instance = 0; \
903  } \
904  TYPE_t **val =(TYPE_t**)(thisobj+fOffset); \
905  return T((val[sub_instance])[index]); \
906  }
907 
908  // pointer to an array of basic types array[n]
923 #if defined(_MSC_VER) && (_MSC_VER <= 1200)
925 #else
927 #endif
928 
929  case kOther_t:
930  default: return 0;
931  }
932 }
933 
934 ////////////////////////////////////////////////////////////////////////////////
935 /// \class TFormLeafInfoDirect
936 /// A small helper class to implement reading a data
937 /// member on an object stored in a TTree.
938 
939 ////////////////////////////////////////////////////////////////////////////////
940 /// Constructor.
941 
943  TFormLeafInfo(from->GetInfo()->GetClass(),0,
944  from->GetInfo()->GetElement(from->GetID()))
945 {
946 }
947 
948 ////////////////////////////////////////////////////////////////////////////////
949 /// Copy this object and its content.
950 
952 {
953  return new TFormLeafInfoDirect(*this);
954 }
955 
956 ////////////////////////////////////////////////////////////////////////////////
957 /// Read the value at the given memory location
958 
959 Double_t TFormLeafInfoDirect::ReadValue(char * /*where*/, Int_t /*instance*/)
960 {
961  Error("ReadValue","Should not be used in a TFormLeafInfoDirect");
962  return 0;
963 }
964 
965 ////////////////////////////////////////////////////////////////////////////////
966 /// Return the underlying value.
967 
968 template <typename T>
969 T TFormLeafInfoDirect::GetValueImpl(TLeaf *leaf, Int_t instance)
970 {
971  return leaf->GetTypedValue<T>(instance);
972 }
973 
975 
976 ////////////////////////////////////////////////////////////////////////////////
977 /// Return the address of the underlying value.
978 
980 {
981  if (leaf->IsA() != TLeafElement::Class()) {
982  return leaf->GetValuePointer();
983  } else {
984  return GetObjectAddress((TLeafElement*)leaf, instance); // instance might be modified
985  }
986 }
987 
988 ////////////////////////////////////////////////////////////////////////////////
989 /// Note this should probably never be executed.
990 
991 void* TFormLeafInfoDirect::GetLocalValuePointer(char *thisobj, Int_t instance)
992 {
993  return TFormLeafInfo::GetLocalValuePointer(thisobj,instance);
994 }
995 
996 ////////////////////////////////////////////////////////////////////////////////
997 /// \class TFormLeafInfoNumerical
998 /// A small helper class to implement reading a numerical value inside a collection
999 
1000 ////////////////////////////////////////////////////////////////////////////////
1001 /// Constructor.
1002 
1004  TFormLeafInfo(0,0,0),
1005  fKind(kind), fIsBool(kFALSE)
1006 {
1007  fElement = new TStreamerElement("data","in collection", 0, fKind, "");
1008 }
1009 
1010 ////////////////////////////////////////////////////////////////////////////////
1011 /// Construct a TFormLeafInfo for the numerical type contained in the collection.
1012 
1014  TFormLeafInfo(0,0,0),
1016 {
1017  if (collection) {
1018  fKind = (EDataType)collection->GetType();
1020  // Could be a bool
1021  if (strcmp( collection->GetCollectionClass()->GetName(), "vector<bool>") == 0
1022  || strncmp( collection->GetCollectionClass()->GetName(), "bitset<", strlen("bitset<") ) ==0 ) {
1023  fIsBool = kTRUE;
1024  fKind = (EDataType)18;
1025  }
1026  }
1027  }
1028  fElement = new TStreamerElement("data","in collection", 0, fKind, "");
1029 }
1030 
1031 ////////////////////////////////////////////////////////////////////////////////
1032 /// Constructor.
1033 
1035  TFormLeafInfo(orig),
1036  fKind(orig.fKind), fIsBool(kFALSE)
1037 {
1038  fElement = new TStreamerElement("data","in collection", 0, fKind, "");
1039 }
1040 
1041 ////////////////////////////////////////////////////////////////////////////////
1042 /// Exception safe swap.
1043 
1045 {
1046  TFormLeafInfo::Swap(other);
1047  std::swap(fKind,other.fKind);
1048  std::swap(fIsBool,other.fIsBool);
1049 }
1050 
1051 ////////////////////////////////////////////////////////////////////////////////
1052 /// Exception safe assignment operator.
1053 
1055 {
1056  TFormLeafInfoNumerical tmp(other);
1057  Swap(tmp);
1058  return *this;
1059 }
1060 
1061 ////////////////////////////////////////////////////////////////////////////////
1062 /// Copy the object and all its content.
1063 
1065 {
1066  return new TFormLeafInfoNumerical(*this);
1067 }
1068 
1069 ////////////////////////////////////////////////////////////////////////////////
1070 /// Destructor
1071 
1073 {
1074  delete fElement;
1075 }
1076 
1077 ////////////////////////////////////////////////////////////////////////////////
1078 /// Return true if the underlying data is a string
1079 
1081 {
1082  if (fIsBool) return kFALSE;
1083  return TFormLeafInfo::IsString();
1084 }
1085 
1086 ////////////////////////////////////////////////////////////////////////////////
1087 /// We reloading all cached information in case the underlying class
1088 /// information has changed (for example when changing from the 'emulated'
1089 /// class to the real class.
1090 
1092 {
1093  //R__ASSERT(fNext==0);
1094 
1095  if (fCounter) return fCounter->Update();
1096  return kFALSE;
1097 }
1098 
1099 namespace {
1100  TStreamerElement *R__GetFakeClonesElem() {
1101  static TStreamerElement gFakeClonesElem("begin","fake",0,
1103  "TClonesArray");
1104  return &gFakeClonesElem;
1105  }
1106 }
1107 
1108 ////////////////////////////////////////////////////////////////////////////////
1109 /// \class TFormLeafInfoClones
1110 /// A small helper class to implement reading a data member
1111 /// on a TClonesArray object stored in a TTree.
1112 
1113 ////////////////////////////////////////////////////////////////////////////////
1114 /// Constructor.
1115 
1117  TFormLeafInfo(classptr,offset,R__GetFakeClonesElem()),fTop(kFALSE)
1118 {
1119 }
1120 
1121 ////////////////////////////////////////////////////////////////////////////////
1122 /// Constructor.
1123 
1125  Bool_t top) :
1126  TFormLeafInfo(classptr,offset,R__GetFakeClonesElem()),fTop(top)
1127 {
1128 }
1129 
1130 ////////////////////////////////////////////////////////////////////////////////
1131 /// Constructor.
1132 
1134  TStreamerElement* element,
1135  Bool_t top) :
1136  TFormLeafInfo(classptr,offset,element),fTop(top)
1137 {
1138 }
1139 
1140 ////////////////////////////////////////////////////////////////////////////////
1141 /// Deep Copy constructor.
1142 
1144  TFormLeafInfo(orig), fTop(orig.fTop)
1145 {
1146 }
1147 
1148 ////////////////////////////////////////////////////////////////////////////////
1149 /// Exception safe swap.
1150 
1152 {
1153  TFormLeafInfo::Swap(other);
1154  std::swap(fTop,other.fTop);
1155 }
1156 
1157 ////////////////////////////////////////////////////////////////////////////////
1158 /// Exception safe assignement operator
1159 
1161 {
1162  TFormLeafInfoClones tmp(orig);
1163  Swap(tmp);
1164  return *this;
1165 }
1166 
1167 ////////////////////////////////////////////////////////////////////////////////
1168 /// Return the current size of the the TClonesArray
1169 
1171 {
1172  if (!fCounter) {
1173  TClass *clonesClass = TClonesArray::Class();
1174  Int_t c_offset = 0;
1175  TStreamerElement *counter = ((TStreamerInfo*)clonesClass->GetStreamerInfo())->GetStreamerElement("fLast",c_offset);
1176  fCounter = new TFormLeafInfo(clonesClass,c_offset,counter);
1177  }
1178  return (Int_t)fCounter->ReadValue((char*)GetLocalValuePointer(leaf)) + 1;
1179 }
1180 
1181 ////////////////////////////////////////////////////////////////////////////////
1182 /// Return the current size of the the TClonesArray
1183 
1185 {
1186  if (!fCounter) {
1187  TClass *clonesClass = TClonesArray::Class();
1188  Int_t c_offset = 0;
1189  TStreamerElement *counter = ((TStreamerInfo*)clonesClass->GetStreamerInfo())->GetStreamerElement("fLast",c_offset);
1190  fCounter = new TFormLeafInfo(clonesClass,c_offset,counter);
1191  }
1192  return (Int_t)fCounter->ReadValue(where) + 1;
1193 }
1194 
1195 ////////////////////////////////////////////////////////////////////////////////
1196 /// Return the value of the underlying data member inside the
1197 /// clones array.
1198 
1199 template <typename T>
1200 T TFormLeafInfoClones::ReadValueImpl(char *where, Int_t instance)
1201 {
1202  if (fNext==0) return 0;
1203  Int_t len,index,sub_instance;
1204  len = fNext->GetArrayLength();
1205  if (len) {
1206  index = instance / len;
1207  sub_instance = instance % len;
1208  } else {
1209  index = instance;
1210  sub_instance = 0;
1211  }
1212  TClonesArray * clones = (TClonesArray*)where;
1213  if (!clones) return 0;
1214  // Note we take advantage of having only one physically variable
1215  // dimension:
1216  char * obj = (char*)clones->UncheckedAt(index);
1217  return fNext->ReadTypedValue<T>(obj,sub_instance);
1218 }
1219 
1222 
1223 ////////////////////////////////////////////////////////////////////////////////
1224 /// Return the pointer to the clonesArray
1225 
1227 {
1228  TClonesArray * clones;
1229  if (fTop) {
1230  if (leaf->InheritsFrom(TLeafObject::Class()) ) {
1231  clones = (TClonesArray*)((TLeafObject*)leaf)->GetObject();
1232  } else {
1233  clones = (TClonesArray*)((TBranchElement*)leaf->GetBranch())->GetObject();
1234  }
1235  } else {
1237  }
1238  return clones;
1239 }
1240 
1241 ////////////////////////////////////////////////////////////////////////////////
1242 /// Return the address of the underlying current value
1243 
1245 {
1246  return TFormLeafInfo::GetLocalValuePointer(where,instance);
1247 }
1248 
1249 ////////////////////////////////////////////////////////////////////////////////
1250 /// Return the value of the underlying data member inside the
1251 /// clones array.
1252 
1253 template <typename T>
1254 T TFormLeafInfoClones::GetValueImpl(TLeaf *leaf, Int_t instance)
1255 {
1256  if (fNext==0) return 0;
1257  Int_t len,index,sub_instance;
1258  len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1259  Int_t primary = fNext->GetPrimaryIndex();
1260  if (len) {
1261  index = instance / len;
1262  sub_instance = instance % len;
1263  } else if (primary>=0) {
1264  index = primary;
1265  sub_instance = instance;
1266  } else {
1267  index = instance;
1268  sub_instance = 0;
1269  }
1271  if (clones==0) return 0;
1272 
1273  // Note we take advantage of having only one physically variable
1274  // dimension:
1275  char * obj = (char*)clones->UncheckedAt(index);
1276  return fNext->ReadTypedValue<T>(obj,sub_instance);
1277 }
1278 
1279 ////////////////////////////////////////////////////////////////////////////////
1280 /// Return the pointer to the clonesArray
1281 
1283 {
1284  TClonesArray * clones = (TClonesArray*)GetLocalValuePointer(leaf);
1285  if (fNext && clones) {
1286  // Same as in TFormLeafInfoClones::GetValue
1287  Int_t len,index,sub_instance;
1288  len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1289  if (len) {
1290  index = instance / len;
1291  sub_instance = instance % len;
1292  } else {
1293  index = instance;
1294  sub_instance = 0;
1295  }
1296  return fNext->GetValuePointer((char*)clones->UncheckedAt(index),
1297  sub_instance);
1298  }
1299  return clones;
1300 }
1301 
1302 ////////////////////////////////////////////////////////////////////////////////
1303 /// Return the pointer to the clonesArray
1304 
1305 void * TFormLeafInfoClones::GetValuePointer(char *where, Int_t instance)
1306 {
1307  TClonesArray * clones = (TClonesArray*) where;
1308  if (fNext) {
1309  // Same as in TFormLeafInfoClones::GetValue
1310  Int_t len,index,sub_instance;
1311  len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1312  if (len) {
1313  index = instance / len;
1314  sub_instance = instance % len;
1315  } else {
1316  index = instance;
1317  sub_instance = 0;
1318  }
1319  return fNext->GetValuePointer((char*)clones->UncheckedAt(index),
1320  sub_instance);
1321  }
1322  return clones;
1323 }
1324 
1325 ////////////////////////////////////////////////////////////////////////////////
1326 /// \class TFormLeafInfoCollectionObject
1327 /// A small helper class to implement reading a data member
1328 /// on a TClonesArray object stored in a TTree.
1329 
1330 ////////////////////////////////////////////////////////////////////////////////
1331 /// Constructor.
1332 
1334  TFormLeafInfo(classptr,0,R__GetFakeClonesElem()),fTop(top)
1335 {
1336 }
1337 
1338 ////////////////////////////////////////////////////////////////////////////////
1339 /// Constructor.
1340 
1342  TFormLeafInfo(orig),fTop(orig.fTop)
1343 {
1344 }
1345 
1346 ////////////////////////////////////////////////////////////////////////////////
1347 /// Exception safe swap.
1348 
1350 {
1351  TFormLeafInfo::Swap(other);
1352  std::swap(fTop,other.fTop);
1353 }
1354 
1355 ////////////////////////////////////////////////////////////////////////////////
1356 /// Exception safe assignement operator
1357 
1359 {
1361  Swap(tmp);
1362  return *this;
1363 }
1364 
1365 ////////////////////////////////////////////////////////////////////////////////
1366 /// Return the current size of the the TClonesArray
1367 
1369 {
1370  return 1;
1371 }
1372 
1373 ////////////////////////////////////////////////////////////////////////////////
1374 /// Return the value of the underlying data member inside the
1375 /// clones array.
1376 
1377 Double_t TFormLeafInfoCollectionObject::ReadValue(char * /* where */, Int_t /* instance */)
1378 {
1379  R__ASSERT(0);
1380  return 0;
1381 }
1382 
1383 ////////////////////////////////////////////////////////////////////////////////
1384 /// Return the pointer to the clonesArray
1385 
1387 {
1388  void* collection;
1389  if (fTop) {
1390  if (leaf->InheritsFrom(TLeafObject::Class()) ) {
1391  collection = ((TLeafObject*)leaf)->GetObject();
1392  } else {
1393  collection = ((TBranchElement*)leaf->GetBranch())->GetObject();
1394  }
1395  } else {
1396  collection = TFormLeafInfo::GetLocalValuePointer(leaf);
1397  }
1398  return collection;
1399 }
1400 
1401 ////////////////////////////////////////////////////////////////////////////////
1402 /// Return the address of the underlying current value
1403 
1405 {
1406  return TFormLeafInfo::GetLocalValuePointer(where,instance);
1407 }
1408 
1409 ////////////////////////////////////////////////////////////////////////////////
1410 /// Return the value of the underlying data member inside the
1411 /// clones array.
1412 
1413 template <typename T>
1414 T TFormLeafInfoCollectionObject::GetValueImpl(TLeaf *leaf, Int_t instance)
1415 {
1416  char * obj = (char*)GetLocalValuePointer(leaf);
1417 
1418  if (fNext==0) return 0;
1419  return fNext->ReadTypedValue<T>(obj,instance);
1420 }
1421 
1423 
1424 ////////////////////////////////////////////////////////////////////////////////
1425 /// Return the pointer to the clonesArray
1426 
1428 {
1429  void *collection = GetLocalValuePointer(leaf);
1430  if (fNext) {
1431  return fNext->GetValuePointer((char*)collection,instance);
1432  }
1433  return collection;
1434 }
1435 
1436 ////////////////////////////////////////////////////////////////////////////////
1437 /// Return the pointer to the clonesArray
1438 
1440 {
1441  if (fNext) {
1442  return fNext->GetValuePointer(where,instance);
1443  }
1444  return where;
1445 }
1446 
1447 ////////////////////////////////////////////////////////////////////////////////
1448 /// \class TFormLeafInfoCollection
1449 /// A small helper class to implement reading a data
1450 /// member on a generic collection object stored in a TTree.
1451 
1452 ////////////////////////////////////////////////////////////////////////////////
1453 /// Cosntructor.
1454 
1456  Long_t offset,
1457  TStreamerElement* element,
1458  Bool_t top) :
1459  TFormLeafInfo(classptr,offset,element),
1460  fTop(top),
1461  fCollClass( 0),
1462  fCollProxy( 0),
1463  fLocalElement( 0)
1464 {
1465  if (element) {
1466  fCollClass = element->GetClass();
1467  } else if (classptr) {
1468  fCollClass = classptr;
1469  }
1470  if (fCollClass
1473 
1476  }
1477 }
1478 
1479 ////////////////////////////////////////////////////////////////////////////////
1480 /// Constructor.
1481 
1483  Long_t offset,
1484  TClass* elementclassptr,
1485  Bool_t top) :
1486  TFormLeafInfo(motherclassptr,offset,
1487  new TStreamerElement("collection","in class",
1488  0,
1489  TStreamerInfo::kAny,
1490  elementclassptr
1491  ? elementclassptr->GetName()
1492  : ( motherclassptr
1493  ? motherclassptr->GetName()
1494  : "Unknwon")
1495  ) ),
1496  fTop(top),
1497  fCollClass( 0),
1498  fCollProxy( 0) ,
1500 {
1501  if (elementclassptr) {
1502  fCollClass = elementclassptr;
1503  } else if (motherclassptr) {
1504  fCollClass = motherclassptr;
1505  }
1506  if (fCollClass
1509  {
1512  }
1513 }
1514 
1515 ////////////////////////////////////////////////////////////////////////////////
1516 /// Constructor.
1517 
1519  TFormLeafInfo(),
1520  fTop(kFALSE),
1521  fCollClass( 0),
1522  fCollProxy( 0),
1523  fLocalElement( 0)
1524 {
1525 }
1526 
1527 ////////////////////////////////////////////////////////////////////////////////
1528 /// Constructor.
1529 
1531  TFormLeafInfo(orig),
1532  fTop( orig.fTop),
1533  fCollClass( orig.fCollClass ),
1535  fCollProxy( orig.fCollProxy ? orig.fCollProxy->Generate() : 0 ),
1536  fLocalElement( 0 ) // humm why not initialize it?
1537 {
1538 }
1539 
1540 ////////////////////////////////////////////////////////////////////////////////
1541 /// Exception safe swap.
1542 
1544 {
1545  TFormLeafInfo::Swap(other);
1546  std::swap(fTop,other.fTop);
1547  std::swap(fClass,other.fClass);
1551 }
1552 
1553 ////////////////////////////////////////////////////////////////////////////////
1554 /// Exception safe assignment operator.
1555 
1557 {
1558  TFormLeafInfoCollection tmp(other);
1559  Swap(tmp);
1560  return *this;
1561 }
1562 
1563 ////////////////////////////////////////////////////////////////////////////////
1564 /// Destructor.
1565 
1567 {
1568  delete fCollProxy;
1569  delete fLocalElement;
1570 }
1571 
1572 ////////////////////////////////////////////////////////////////////////////////
1573 /// Copy of the object and its content.
1574 
1576 {
1577  return new TFormLeafInfoCollection(*this);
1578 }
1579 
1580 ////////////////////////////////////////////////////////////////////////////////
1581 /// We reloading all cached information in case the underlying class
1582 /// information has changed (for example when changing from the 'emulated'
1583 /// class to the real class.
1584 
1586 {
1587  Bool_t changed = kFALSE;
1588  TClass * new_class = TClass::GetClass(fCollClassName);
1589  if (new_class!=fCollClass) {
1590  delete fCollProxy; fCollProxy = 0;
1591  fCollClass = new_class;
1594  }
1595  changed = kTRUE;
1596  }
1597  return changed || TFormLeafInfo::Update();
1598 }
1599 
1600 ////////////////////////////////////////////////////////////////////////////////
1601 /// Return true if the underlying data has a array size counter
1602 
1604 {
1605  return fCounter!=0 || fCollProxy!=0;
1606 }
1607 
1608 ////////////////////////////////////////////////////////////////////////////////
1609 /// Return the current size of the the TClonesArray
1610 
1612 {
1613  void *ptr = GetLocalValuePointer(leaf);
1614 
1615  if (fCounter) { return (Int_t)fCounter->ReadValue((char*)ptr); }
1616 
1618  if (ptr==0) return 0;
1620  return (Int_t)fCollProxy->Size();
1621 }
1622 
1623 ////////////////////////////////////////////////////////////////////////////////
1624 /// Return the size of the underlying array for the current entry in the TTree.
1625 
1627 {
1628  if (fCounter) { return (Int_t)fCounter->ReadValue(where); }
1630  if (where==0) return 0;
1631  void *ptr = GetLocalValuePointer(where,0);
1633  return (Int_t)fCollProxy->Size();
1634 }
1635 
1636 ////////////////////////////////////////////////////////////////////////////////
1637 /// Return the current size of the the TClonesArray
1638 
1640 {
1641  void *ptr = GetLocalValuePointer(leaf,instance);
1642  if (fCounter) {
1643  return (Int_t)fCounter->ReadValue((char*)ptr);
1644  }
1646  if (ptr==0) return 0;
1648  return (Int_t)fCollProxy->Size();
1649 }
1650 
1651 ////////////////////////////////////////////////////////////////////////////////
1652 /// Return the value of the underlying data member inside the
1653 /// clones array.
1654 
1655 template <typename T>
1656 T TFormLeafInfoCollection::ReadValueImpl(char *where, Int_t instance)
1657 {
1658  if (fNext==0) return 0;
1659  UInt_t len,index,sub_instance;
1660  len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1661  Int_t primary = fNext->GetPrimaryIndex();
1662  if (len) {
1663  index = instance / len;
1664  sub_instance = instance % len;
1665  } else if (primary>=0) {
1666  index = primary;
1667  sub_instance = instance;
1668  } else {
1669  index = instance;
1670  sub_instance = 0;
1671  }
1672 
1674  void *ptr = GetLocalValuePointer(where,instance);
1676 
1677  // Note we take advantage of having only one physically variable
1678  // dimension:
1679 
1680  char * obj = (char*)fCollProxy->At(index);
1681  if (fCollProxy->HasPointers()) obj = *(char**)obj;
1682  return fNext->ReadTypedValue<T>(obj,sub_instance);
1683 }
1684 
1687 
1688 ////////////////////////////////////////////////////////////////////////////////
1689 /// Return the pointer to the clonesArray
1690 
1692 {
1693  void *collection;
1694  if (fTop) {
1695  if (leaf->InheritsFrom(TLeafObject::Class()) ) {
1696  collection = ((TLeafObject*)leaf)->GetObject();
1697  } else {
1698  collection = ((TBranchElement*)leaf->GetBranch())->GetObject();
1699  }
1700  } else {
1701  collection = TFormLeafInfo::GetLocalValuePointer(leaf);
1702  }
1703  return collection;
1704 }
1705 
1706 ////////////////////////////////////////////////////////////////////////////////
1707 /// Return the address of the local value
1708 
1710 {
1711  return TFormLeafInfo::GetLocalValuePointer(where,instance);
1712 }
1713 
1714 ////////////////////////////////////////////////////////////////////////////////
1715 /// Return the value of the underlying data member inside the
1716 /// clones array.
1717 
1718 template <typename T>
1719 T TFormLeafInfoCollection::GetValueImpl(TLeaf *leaf, Int_t instance)
1720 {
1721  if (fNext==0) return 0;
1722  Int_t len,index,sub_instance;
1723  len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1724  Int_t primary = fNext->GetPrimaryIndex();
1725  if (len) {
1726  index = instance / len;
1727  sub_instance = instance % len;
1728  } else if (primary>=0) {
1729  index = primary;
1730  sub_instance = instance;
1731  } else {
1732  index = instance;
1733  sub_instance = 0;
1734  }
1735 
1737  void *coll = GetLocalValuePointer(leaf);
1739 
1740  // Note we take advantage of having only one physically variable
1741  // dimension:
1742  char * obj = (char*)fCollProxy->At(index);
1743  if (obj==0) return 0;
1744  if (fCollProxy->HasPointers()) obj = *(char**)obj;
1745  if (obj==0) return 0;
1746  return fNext->ReadTypedValue<T>(obj,sub_instance);
1747 }
1748 
1749 ////////////////////////////////////////////////////////////////////////////////
1750 /// Return the pointer to the clonesArray
1751 
1753 {
1755 
1756  void *collection = GetLocalValuePointer(leaf);
1757 
1758  if (fNext) {
1759  // Same as in TFormLeafInfoClones::GetValue
1760  Int_t len,index,sub_instance;
1761  if (fNext->fElement &&
1762  (fNext->fNext || !fNext->IsString()) ) {
1763  len = fNext->GetArrayLength();
1764  } else {
1765  len = 0;
1766  }
1767  if (len) {
1768  index = instance / len;
1769  sub_instance = instance % len;
1770  } else {
1771  index = instance;
1772  sub_instance = 0;
1773  }
1774  TVirtualCollectionProxy::TPushPop helper(fCollProxy,collection);
1775  char * obj = (char*)fCollProxy->At(index);
1776  if (fCollProxy->HasPointers()) obj = *(char**)obj;
1777  return fNext->GetValuePointer(obj,sub_instance);
1778  }
1779  return collection;
1780 }
1781 
1782 ////////////////////////////////////////////////////////////////////////////////
1783 /// Return the pointer to the clonesArray
1784 
1786 {
1788 
1789  void *collection = where;
1790 
1791  if (fNext) {
1792  // Same as in TFormLeafInfoClones::GetValue
1793  Int_t len,index,sub_instance;
1794  len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1795  if (len) {
1796  index = instance / len;
1797  sub_instance = instance % len;
1798  } else {
1799  index = instance;
1800  sub_instance = 0;
1801  }
1802  TVirtualCollectionProxy::TPushPop helper(fCollProxy,collection);
1803  char * obj = (char*)fCollProxy->At(index);
1804  if (fCollProxy->HasPointers()) obj = *(char**)obj;
1805  return fNext->GetValuePointer(obj,sub_instance);
1806  }
1807  return collection;
1808 }
1809 
1810 ////////////////////////////////////////////////////////////////////////////////
1811 /// \class TFormLeafInfoCollectionSize
1812 /// Used to return the size of a collection
1813 
1814 ////////////////////////////////////////////////////////////////////////////////
1815 /// Constructor.
1816 
1818  TFormLeafInfo(), fCollClass(classptr), fCollProxy(0)
1819 {
1820  if (fCollClass
1823 
1826  }
1827 }
1828 
1829 ////////////////////////////////////////////////////////////////////////////////
1830 /// Constructor.
1831 
1833  TClass* classptr,Long_t offset,TStreamerElement* element) :
1834  TFormLeafInfo(classptr,offset,element), fCollClass(element->GetClassPointer()), fCollProxy(0)
1835 {
1836  if (fCollClass
1839 
1842  }
1843 }
1844 
1845 ////////////////////////////////////////////////////////////////////////////////
1846 /// Constructor.
1847 
1850 {
1851 }
1852 
1853 ////////////////////////////////////////////////////////////////////////////////
1854 /// Constructor.
1855 
1857  const TFormLeafInfoCollectionSize& orig) : TFormLeafInfo(),
1858  fCollClass(orig.fCollClass),
1860  fCollProxy(orig.fCollProxy?orig.fCollProxy->Generate():0)
1861 {
1862 }
1863 
1864 ////////////////////////////////////////////////////////////////////////////////
1865 /// Exception safe assignment operator.
1866 
1868 {
1869  TFormLeafInfoCollectionSize tmp(other);
1870  Swap(tmp);
1871  return *this;
1872 }
1873 
1874 ////////////////////////////////////////////////////////////////////////////////
1875 /// Exception safe swap.
1876 
1878 {
1879  TFormLeafInfo::Swap(other);
1883 }
1884 
1885 ////////////////////////////////////////////////////////////////////////////////
1886 /// Destructor.
1887 
1889 {
1890  delete fCollProxy;
1891 }
1892 
1893 ////////////////////////////////////////////////////////////////////////////////
1894 /// Copy the object and all of its content.
1895 
1897 {
1898  return new TFormLeafInfoCollectionSize(*this);
1899 }
1900 
1901 ////////////////////////////////////////////////////////////////////////////////
1902 /// We reloading all cached information in case the underlying class
1903 /// information has changed (for example when changing from the 'emulated'
1904 /// class to the real class.
1905 
1907 {
1908  Bool_t changed = kFALSE;
1909  TClass *new_class = TClass::GetClass(fCollClassName);
1910  if (new_class!=fCollClass) {
1911  delete fCollProxy; fCollProxy = 0;
1912  fCollClass = new_class;
1915  }
1916  changed = kTRUE;
1917  }
1918  return changed;
1919 }
1920 
1921 ////////////////////////////////////////////////////////////////////////////////
1922 /// Not implemented.
1923 
1924 void *TFormLeafInfoCollectionSize::GetValuePointer(TLeaf * /* leaf */, Int_t /* instance */)
1925 {
1926  Error("GetValuePointer","This should never be called");
1927  return 0;
1928 }
1929 
1930 ////////////////////////////////////////////////////////////////////////////////
1931 /// Not implemented.
1932 
1933 void *TFormLeafInfoCollectionSize::GetValuePointer(char * /* from */, Int_t /* instance */)
1934 {
1935  Error("GetValuePointer","This should never be called");
1936  return 0;
1937 }
1938 
1939 ////////////////////////////////////////////////////////////////////////////////
1940 /// Not implemented.
1941 
1943 {
1944  Error("GetLocalValuePointer","This should never be called");
1945  return 0;
1946 }
1947 
1948 ////////////////////////////////////////////////////////////////////////////////
1949 /// Not implemented.
1950 
1951 void *TFormLeafInfoCollectionSize::GetLocalValuePointer( char * /* from */, Int_t /* instance */)
1952 {
1953  Error("GetLocalValuePointer","This should never be called");
1954  return 0;
1955 }
1956 
1957 ////////////////////////////////////////////////////////////////////////////////
1958 /// Return the value of the underlying pointer data member
1959 
1961 {
1963  if (where==0) return 0;
1964  void *ptr = fElement ? TFormLeafInfo::GetLocalValuePointer(where) : where;
1966  return (Int_t)fCollProxy->Size();
1967 }
1968 
1969 ////////////////////////////////////////////////////////////////////////////////
1970 /// \class TFormLeafInfoPointer
1971 /// A small helper class to implement reading a data
1972 /// member by following a pointer inside a branch of TTree.
1973 
1974 ////////////////////////////////////////////////////////////////////////////////
1975 /// Constructor.
1976 
1978  Long_t offset,
1979  TStreamerElement* element) :
1980  TFormLeafInfo(classptr,offset,element)
1981 {
1982 }
1983 
1984 ////////////////////////////////////////////////////////////////////////////////
1985 /// Copy the object and all of its contnet.
1986 
1988 {
1989  return new TFormLeafInfoPointer(*this);
1990 }
1991 
1992 ////////////////////////////////////////////////////////////////////////////////
1993 /// Return the value of the underlying pointer data member
1994 
1995 template <typename T>
1996 T TFormLeafInfoPointer::ReadValueImpl(char *where, Int_t instance)
1997 {
1998  if (!fNext) return 0;
1999  char * whereoffset = where+fOffset;
2000  switch (fElement->GetNewType()) {
2001  // basic types
2004  case TStreamerInfo::kAnyp:
2005  case TStreamerInfo::kAnyP:
2006  case TStreamerInfo::kSTLp:
2007  {TObject **obj = (TObject**)(whereoffset);
2008  return obj && *obj ? fNext->ReadTypedValue<T>((char*)*obj,instance) : 0; }
2009 
2014  case TStreamerInfo::kAny:
2015  case TStreamerInfo::kBase:
2016  case TStreamerInfo::kSTL:
2017  {
2018  TObject *obj = (TObject*)(whereoffset);
2019  return fNext->ReadTypedValue<T>((char*)obj,instance);
2020  }
2021 
2025  {
2026  Int_t len, index, sub_instance;
2027 
2028  if (fNext) len = fNext->GetArrayLength();
2029  else len = 1;
2030  if (len) {
2031  index = instance / len;
2032  sub_instance = instance % len;
2033  } else {
2034  index = instance;
2035  sub_instance = 0;
2036  }
2037 
2038  whereoffset += index*fElement->GetClassPointer()->Size();
2039 
2040  TObject *obj = (TObject*)(whereoffset);
2041  return fNext->ReadTypedValue<T>((char*)obj,sub_instance);
2042  }
2043 
2049  {
2050  TObject *obj = (TObject*)(whereoffset);
2051  return fNext->ReadTypedValue<T>((char*)obj,instance);
2052  }
2053 
2054  case kOther_t:
2055  default: return 0;
2056  }
2057 
2058 }
2059 
2060 ////////////////////////////////////////////////////////////////////////////////
2061 /// Return the value of the underlying pointer data member
2062 
2063 template <typename T>
2064 T TFormLeafInfoPointer::GetValueImpl(TLeaf *leaf, Int_t instance)
2065 {
2066  if (!fNext) return 0;
2067  char * where = (char*)GetLocalValuePointer(leaf,instance);
2068  if (where==0) return 0;
2069  return fNext->ReadTypedValue<T>(where,instance);
2070 }
2071 
2074 
2075 ////////////////////////////////////////////////////////////////////////////////
2076 /// \class TFormLeafInfoMethod
2077 /// Asmall helper class to implement executing a method
2078 /// of an object stored in a TTree
2079 
2080 ////////////////////////////////////////////////////////////////////////////////
2081 /// Constructor.
2082 
2084  TMethodCall *method) :
2085  TFormLeafInfo(classptr,0,0),fMethod(method),
2086  fResult(0), fCopyFormat(),fDeleteFormat(),fValuePointer(0),fIsByValue(kFALSE)
2087 {
2088  if (method) {
2089  fMethodName = method->GetMethodName();
2090  fParams = method->GetParams();
2092  if (r == TMethodCall::kOther) {
2093  const char* rtype = fMethod->GetMethod()->GetReturnTypeName();
2094  Long_t rprop = fMethod->GetMethod()->Property();
2095  if (rtype[strlen(rtype)-1]!='*' &&
2096  rtype[strlen(rtype)-1]!='&' &&
2097  !(rprop & (kIsPointer|kIsReference)) ) {
2098  fCopyFormat = "new ";
2099  fCopyFormat += rtype;
2100  fCopyFormat += "(*(";
2101  fCopyFormat += rtype;
2102  fCopyFormat += "*)0x%lx)";
2103 
2104  fDeleteFormat = "delete (";
2105  fDeleteFormat += rtype;
2106  fDeleteFormat += "*)0x%lx";
2107 
2108  fIsByValue = kTRUE;
2109  }
2110  }
2111  }
2112 }
2113 
2114 ////////////////////////////////////////////////////////////////////////////////
2115 /// Constructor.
2116 
2118  : TFormLeafInfo(orig)
2119 {
2120  fMethodName = orig.fMethodName;
2121  fParams = orig.fParams ;
2122  fResult = orig.fResult;
2123  if (orig.fMethod) {
2124  fMethod = new TMethodCall();
2125  fMethod->Init(orig.fMethod->GetMethod());
2126  } else {
2127  fMethod = 0;
2128  }
2129  fCopyFormat = orig.fCopyFormat;
2131  fValuePointer = 0;
2132  fIsByValue = orig.fIsByValue;
2133 }
2134 
2135 ////////////////////////////////////////////////////////////////////////////////
2136 /// Exception safe swap.
2137 
2139 {
2140  TFormLeafInfo::Swap(other);
2141  std::swap(fMethod,other.fMethod);
2143  std::swap(fParams,other.fParams);
2144  std::swap(fResult,other.fResult);
2149 }
2150 
2151 ////////////////////////////////////////////////////////////////////////////////
2152 /// Exception safe assignment operator.
2153 
2155 {
2156  TFormLeafInfoMethod tmp(other);
2157  Swap(tmp);
2158  return *this;
2159 }
2160 
2161 ////////////////////////////////////////////////////////////////////////////////
2162 /// Destructor.
2163 
2165 {
2166  if (fValuePointer) {
2168  }
2169  delete fMethod;
2170 }
2171 
2172 ////////////////////////////////////////////////////////////////////////////////
2173 /// Copy the object and all its content.
2174 
2176 {
2177  return new TFormLeafInfoMethod(*this);
2178 }
2179 
2180 ////////////////////////////////////////////////////////////////////////////////
2181 /// Return the type of the underlying return value
2182 
2184 {
2185  if (fNext) return fNext->GetClass();
2187  if (r!=TMethodCall::kOther) return 0;
2188  TString return_type = gInterpreter->TypeName(fMethod->GetMethod()->GetReturnTypeName());
2189  return TClass::GetClass(return_type.Data());
2190 }
2191 
2192 ////////////////////////////////////////////////////////////////////////////////
2193 /// Return true if the return value is integral.
2194 
2196 {
2198  if (r == TMethodCall::kLong) {
2199  return kTRUE;
2200  } else return kFALSE;
2201 }
2202 
2203 ////////////////////////////////////////////////////////////////////////////////
2204 /// Return true if the return value is a string.
2205 
2207 {
2208  if (fNext) return fNext->IsString();
2209 
2211  return (r==TMethodCall::kString);
2212 }
2213 
2214 ////////////////////////////////////////////////////////////////////////////////
2215 /// We reloading all cached information in case the underlying class
2216 /// information has changed (for example when changing from the 'emulated'
2217 /// class to the real class.
2218 
2220 {
2221  if (!TFormLeafInfo::Update()) return kFALSE;
2222  delete fMethod;
2224  return kTRUE;
2225 }
2226 
2227 ////////////////////////////////////////////////////////////////////////////////
2228 /// This is implemented here because some compiler want ALL the
2229 /// signature of an overloaded function to be re-implemented.
2230 
2232  Int_t instance)
2233 {
2234  return TFormLeafInfo::GetLocalValuePointer( from, instance);
2235 }
2236 
2237 ////////////////////////////////////////////////////////////////////////////////
2238 /// Return the address of the lcoal underlying value.
2239 
2241  Int_t /*instance*/)
2242 {
2243  void *thisobj = from;
2244  if (!thisobj) return 0;
2245 
2247  fResult = 0;
2248 
2249  if (r == TMethodCall::kLong) {
2250  Long_t l = 0;
2251  fMethod->Execute(thisobj, l);
2252  fResult = (Double_t) l;
2253  // Get rid of temporary return object.
2254  gInterpreter->ClearStack();
2255  return &fResult;
2256 
2257  } else if (r == TMethodCall::kDouble) {
2258  Double_t d = 0;
2259  fMethod->Execute(thisobj, d);
2260  fResult = (Double_t) d;
2261  // Get rid of temporary return object.
2262  gInterpreter->ClearStack();
2263  return &fResult;
2264 
2265  } else if (r == TMethodCall::kString) {
2266  char *returntext = 0;
2267  fMethod->Execute(thisobj,&returntext);
2268  gInterpreter->ClearStack();
2269  return returntext;
2270 
2271  } else if (r == TMethodCall::kOther) {
2272  char * char_result = 0;
2273  if (fIsByValue) {
2274  if (fValuePointer) {
2275  gROOT->ProcessLine(Form(fDeleteFormat.Data(),fValuePointer));
2276  fValuePointer = 0;
2277  }
2278  }
2279  fMethod->Execute(thisobj, &char_result);
2280  if (fIsByValue) {
2281  fValuePointer = (char*)gInterpreter->Calc(Form(fCopyFormat.Data(),char_result));
2282  char_result = (char*)fValuePointer;
2283  }
2284  gInterpreter->ClearStack();
2285  return char_result;
2286 
2287  }
2288  return 0;
2289 }
2290 
2291 ////////////////////////////////////////////////////////////////////////////////
2292 /// Execute the method on the given address
2293 
2294 template <typename T>
2295 T TFormLeafInfoMethod::ReadValueImpl(char *where, Int_t instance)
2296 {
2297  void *thisobj = where;
2298  if (!thisobj) return 0;
2299 
2301  T result = 0;
2302 
2303  if (r == TMethodCall::kLong) {
2304  Long_t l = 0;
2305  fMethod->Execute(thisobj, l);
2306  result = (T) l;
2307 
2308  } else if (r == TMethodCall::kDouble) {
2309  Double_t d = 0;
2310  fMethod->Execute(thisobj, d);
2311  result = (T) d;
2312 
2313  } else if (r == TMethodCall::kString) {
2314  char *returntext = 0;
2315  fMethod->Execute(thisobj,&returntext);
2316  result = T((Long_t) returntext);
2317 
2318  } else if (fNext) {
2319  char * char_result = 0;
2320  fMethod->Execute(thisobj, &char_result);
2321  result = fNext->ReadTypedValue<T>(char_result,instance);
2322 
2323  } else fMethod->Execute(thisobj);
2324 
2325  // Get rid of temporary return object.
2326  gInterpreter->ClearStack();
2327  return result;
2328 }
2329 
2331 
2332 ////////////////////////////////////////////////////////////////////////////////
2333 /// \class TFormLeafInfoMultiVarDim
2334 /// A helper class to implement reading a
2335 /// data member on a variable size array inside a TClonesArray object stored in
2336 /// a TTree. This is the version used when the data member is inside a
2337 /// non-split object.
2338 
2339 ////////////////////////////////////////////////////////////////////////////////
2340 /// Constructor.
2341 
2343  Long_t offset,
2344  TStreamerElement* element,
2345  TFormLeafInfo* parent) :
2346  TFormLeafInfo(classptr,offset,element),fNsize(0),fCounter2(0),fSumOfSizes(0),
2347  fDim(0),fVirtDim(-1),fPrimaryIndex(-1),fSecondaryIndex(-1)
2348 {
2349  if (element && element->InheritsFrom(TStreamerBasicPointer::Class())) {
2350  TStreamerBasicPointer * elem = (TStreamerBasicPointer*)element;
2351 
2352  Int_t counterOffset = 0;
2353  TStreamerElement* counter = ((TStreamerInfo*)classptr->GetStreamerInfo())->GetStreamerElement(elem->GetCountName(),counterOffset);
2354  if (!parent) return;
2355  fCounter2 = parent->DeepCopy();
2356  TFormLeafInfo ** next = &(fCounter2->fNext);
2357  while(*next != 0) next = &( (*next)->fNext);
2358  *next = new TFormLeafInfo(classptr,counterOffset,counter);
2359 
2360  } else Error("Constructor","Called without a proper TStreamerElement");
2361 }
2362 
2363 ////////////////////////////////////////////////////////////////////////////////
2364 /// Constructor.
2365 
2367  TFormLeafInfo(0,0,0),fNsize(0),fCounter2(0),fSumOfSizes(0),
2369 {
2370 }
2371 
2372 ////////////////////////////////////////////////////////////////////////////////
2373 /// Constructor.
2374 
2376 {
2377  fNsize = orig.fNsize;
2378  fSizes.Copy(fSizes);
2379  fCounter2 = orig.fCounter2?orig.fCounter2->DeepCopy():0;
2380  fSumOfSizes = orig.fSumOfSizes;
2381  fDim = orig.fDim;
2382  fVirtDim = orig.fVirtDim;
2385 }
2386 
2387 ////////////////////////////////////////////////////////////////////////////////
2388 /// Exception safe swap.
2389 
2391 {
2392  TFormLeafInfo::Swap(other);
2393  std::swap(fNsize,other.fNsize);
2394  std::swap(fSizes,other.fSizes);
2396  std::swap(fDim,other.fDim);
2397  std::swap(fVirtDim,other.fVirtDim);
2400 }
2401 
2402 ////////////////////////////////////////////////////////////////////////////////
2403 /// Exception safe assignment operator.
2404 
2406 {
2407  TFormLeafInfoMultiVarDim tmp(other);
2408  Swap(tmp);
2409  return *this;
2410 }
2411 
2412 ////////////////////////////////////////////////////////////////////////////////
2413 /// Copy the object and all its content.
2414 
2416 {
2417  return new TFormLeafInfoMultiVarDim(*this);
2418 }
2419 
2420 ////////////////////////////////////////////////////////////////////////////////
2421 /// Destructor.
2422 
2424 {
2425  delete fCounter2;
2426 }
2427 
2428 /* The proper indexing and unwinding of index is done by prior leafinfo in the chain. */
2429 //virtual Double_t TFormLeafInfoMultiVarDim::ReadValue(char *where, Int_t instance = 0) {
2430 // return TFormLeafInfo::ReadValue(where,instance);
2431 //}
2432 
2433 ////////////////////////////////////////////////////////////////////////////////
2434 /// Load the current array sizes.
2435 
2437 {
2438  if (fElement) {
2439  TLeaf *leaf = (TLeaf*)branch->GetListOfLeaves()->At(0);
2440  if (fCounter) fNsize = (Int_t)fCounter->GetValue(leaf);
2441  else fNsize = fCounter2->GetCounterValue(leaf);
2442  if (fNsize > fSizes.GetSize()) fSizes.Set(fNsize);
2443  fSumOfSizes = 0;
2444  for (Int_t i=0; i<fNsize; i++) {
2445  Int_t size = (Int_t)fCounter2->GetValue(leaf,i);
2446  fSumOfSizes += size;
2447  fSizes.AddAt( size, i );
2448  }
2449  return;
2450  }
2451  if (!fCounter2 || !fCounter) return;
2452  TBranchElement *br = dynamic_cast<TBranchElement*>(branch);
2453  R__ASSERT(br);
2454  fNsize = br->GetBranchCount()->GetNdata();
2455  if (fNsize > fSizes.GetSize()) fSizes.Set(fNsize);
2456  fSumOfSizes = 0;
2457  for (Int_t i=0; i<fNsize; i++) {
2458  Int_t size = (Int_t)fCounter2->GetValue((TLeaf*)br->GetBranchCount2()->GetListOfLeaves()->At(0),i);
2459  fSumOfSizes += size;
2460  fSizes.AddAt( size, i );
2461  }
2462 }
2463 
2464 ////////////////////////////////////////////////////////////////////////////////
2465 /// Return the index vlaue of the primary index.
2466 
2468 {
2469  return fPrimaryIndex;
2470 }
2471 
2472 ////////////////////////////////////////////////////////////////////////////////
2473 /// Return the size of the requested sub-array.
2474 
2476 {
2477  if (index >= fSizes.GetSize()) {
2478  return -1;
2479  } else {
2480  return fSizes.At(index);
2481  }
2482 }
2483 
2484 ////////////////////////////////////////////////////////////////////////////////
2485 /// Set the current value of the primary index.
2486 
2488 {
2489  fPrimaryIndex = index;
2490 }
2491 
2492 ////////////////////////////////////////////////////////////////////////////////
2493 /// Set the current value of the primary index.
2494 
2496 {
2497  fSecondaryIndex = index;
2498 }
2499 
2500 ////////////////////////////////////////////////////////////////////////////////
2501 /// Set the sizes of the sub-array.
2502 
2504 {
2505  fSumOfSizes += (val - fSizes.At(index));
2506  fSizes.AddAt(val,index);
2507 }
2508 
2509 ////////////////////////////////////////////////////////////////////////////////
2510 /// Get the total size.
2511 
2513 {
2514  return fSumOfSizes;
2515 }
2516 
2517 ////////////////////////////////////////////////////////////////////////////////
2518 
2520  Int_t /*instance*/)
2521 {
2522  /* The proper indexing and unwinding of index need to be done by prior leafinfo in the chain. */
2523  Error("GetValue","This should never be called");
2524  return 0;
2525 }
2526 
2527 
2528 ////////////////////////////////////////////////////////////////////////////////
2529 /// Return the index of the dimension which varies
2530 /// for each elements of an enclosing array (typically a TClonesArray)
2531 
2533 {
2534  return fDim;
2535 }
2536 
2537 ////////////////////////////////////////////////////////////////////////////////
2538 /// Return the virtual index (for this expression) of the dimension which varies
2539 /// for each elements of an enclosing array (typically a TClonesArray)
2540 
2542 {
2543  return fVirtDim;
2544 }
2545 
2546 ////////////////////////////////////////////////////////////////////////////////
2547 /// We reloading all cached information in case the underlying class
2548 /// information has changed (for example when changing from the 'emulated'
2549 /// class to the real class.
2550 
2552 {
2553  Bool_t res = TFormLeafInfo::Update();
2554  if (fCounter2) fCounter2->Update();
2555  return res;
2556 }
2557 
2558 ////////////////////////////////////////////////////////////////////////////////
2559 /// Update the sizes of the arrays.
2560 
2562 {
2563  if (!garr) return;
2564  if (garr->GetSize()<fNsize) garr->Set(fNsize);
2565  for (Int_t i=0; i<fNsize; i++) {
2566  Int_t local = fSizes.At(i);
2567  Int_t global = garr->At(i);
2568  if (global==0 || local<global) global = local;
2569  garr->AddAt(global,i);
2570  }
2571 }
2572 
2573 ////////////////////////////////////////////////////////////////////////////////
2574 /// \class TFormLeafInfoMultiVarDimDirect
2575 /// A small helper class to implement reading
2576 /// a data member on a variable size array inside a TClonesArray object stored
2577 /// in a TTree. This is the version used for split access
2578 
2579 ////////////////////////////////////////////////////////////////////////////////
2580 /// Copy the object and all its content.
2581 
2583 {
2584  return new TFormLeafInfoMultiVarDimDirect(*this);
2585 }
2586 
2587 ////////////////////////////////////////////////////////////////////////////////
2588 /// Return the undersying value.
2589 
2590 template <typename T>
2591 T TFormLeafInfoMultiVarDimDirect::GetValueImpl(TLeaf *leaf, Int_t instance)
2592 {
2593  return ((TLeafElement*)leaf)->GetTypedValueSubArray<T>(fPrimaryIndex,instance);
2594 }
2595 
2597 
2598 ////////////////////////////////////////////////////////////////////////////////
2599 /// Not implemented.
2600 
2602 {
2603  Error("ReadValue","This should never be called");
2604  return 0;
2605 }
2606 
2607 ////////////////////////////////////////////////////////////////////////////////
2608 /// \class TFormLeafInfoMultiVarDimCollection
2609 /// A small helper class to implement reading
2610 /// a data member on a variable size array inside a TClonesArray object stored
2611 /// in a TTree. This is the version used for split access
2612 
2613 ////////////////////////////////////////////////////////////////////////////////
2614 /// Constructor.
2615 
2617  TClass* motherclassptr,
2618  Long_t offset,
2619  TClass* elementclassptr,
2620  TFormLeafInfo *parent) :
2621  TFormLeafInfoMultiVarDim(motherclassptr,offset,
2622  new TStreamerElement("collection","in class",
2623  0,
2624  TStreamerInfo::kAny,
2625  elementclassptr
2626  ? elementclassptr->GetName()
2627  : ( motherclassptr
2628  ? motherclassptr->GetName()
2629  : "Unknwon")
2630  )
2631  )
2632 {
2633  R__ASSERT(parent);
2634  fCounter = parent->DeepCopy();
2635  fCounter2 = parent->DeepCopy();
2636  TFormLeafInfo ** next = &(fCounter2->fNext);
2637  while(*next != 0) next = &( (*next)->fNext);
2638  *next = new TFormLeafInfoCollectionSize(elementclassptr);
2639 }
2640 
2641 ////////////////////////////////////////////////////////////////////////////////
2642 /// Constructor.
2643 
2645  TClass* motherclassptr,
2646  Long_t offset,
2647  TStreamerElement* element,
2648  TFormLeafInfo *parent) :
2649  TFormLeafInfoMultiVarDim(motherclassptr,offset,element)
2650 {
2651  R__ASSERT(parent && element);
2652  fCounter = parent->DeepCopy();
2653  fCounter2 = parent->DeepCopy();
2654  TFormLeafInfo ** next = &(fCounter2->fNext);
2655  while(*next != 0) next = &( (*next)->fNext);
2656  *next = new TFormLeafInfoCollectionSize(motherclassptr,offset,element);
2657 }
2658 
2659 ////////////////////////////////////////////////////////////////////////////////
2660 /// Copy the object and all its content.
2661 
2663 {
2664  return new TFormLeafInfoMultiVarDimCollection(*this);
2665 }
2666 
2667 ////////////////////////////////////////////////////////////////////////////////
2668 
2670  Int_t /* instance */)
2671 {
2672  /* The proper indexing and unwinding of index need to be done by prior leafinfo in the chain. */
2673  Error("GetValue","This should never be called");
2674  return 0;
2675 }
2676 
2677 ////////////////////////////////////////////////////////////////////////////////
2678 /// Load the current array sizes.
2679 
2681 {
2683 
2684  TLeaf *leaf = (TLeaf*)branch->GetListOfLeaves()->At(0);
2686 
2687  if (fNsize > fSizes.GetSize()) fSizes.Set(fNsize);
2688  fSumOfSizes = 0;
2689  for (Int_t i=0; i<fNsize; i++) {
2690  Int_t size = (Int_t)fCounter2->GetValue(leaf,i);
2691  fSumOfSizes += size;
2692  fSizes.AddAt( size, i );
2693  }
2694  return;
2695 }
2696 
2697 ////////////////////////////////////////////////////////////////////////////////
2698 /// Return the value of the underlying data.
2699 
2700 template <typename T>
2701 T TFormLeafInfoMultiVarDimCollection::ReadValueImpl(char *where, Int_t instance)
2702 {
2703  if (fSecondaryIndex>=0) {
2704  UInt_t len = fNext->GetArrayLength();
2705  if (len) {
2706  instance = fSecondaryIndex*len;
2707  } else {
2708  instance = fSecondaryIndex;
2709  }
2710  }
2711  return fNext->ReadTypedValue<T>(where,instance);
2712 }
2713 
2715 
2716 
2717 ////////////////////////////////////////////////////////////////////////////////
2718 /// \class TFormLeafInfoMultiVarDimClones
2719 /// A small helper class to implement reading
2720 /// a data member on a variable size array inside a TClonesArray object stored
2721 /// in a TTree. This is the version used for split access
2722 
2723 ////////////////////////////////////////////////////////////////////////////////
2724 /// Constructor.
2725 
2727  TClass* motherclassptr,
2728  Long_t offset,
2729  TClass* elementclassptr,
2730  TFormLeafInfo *parent) :
2731  TFormLeafInfoMultiVarDim(motherclassptr,offset,
2732  new TStreamerElement("clones","in class",
2733  0,
2734  TStreamerInfo::kAny,
2735  elementclassptr
2736  ? elementclassptr->GetName()
2737  : ( motherclassptr
2738  ? motherclassptr->GetName()
2739  : "Unknwon")
2740  )
2741  )
2742 {
2743  R__ASSERT(parent);
2744  fCounter = parent->DeepCopy();
2745  fCounter2 = parent->DeepCopy();
2746  TFormLeafInfo ** next = &(fCounter2->fNext);
2747  while(*next != 0) next = &( (*next)->fNext);
2748  *next = new TFormLeafInfoClones(elementclassptr);
2749 }
2750 
2751 ////////////////////////////////////////////////////////////////////////////////
2752 /// Constructor.
2753 
2755  TClass* motherclassptr,
2756  Long_t offset,
2757  TStreamerElement* element,
2758  TFormLeafInfo *parent) :
2759  TFormLeafInfoMultiVarDim(motherclassptr,offset,element)
2760 {
2761  R__ASSERT(parent && element);
2762  fCounter = parent->DeepCopy();
2763  fCounter2 = parent->DeepCopy();
2764  TFormLeafInfo ** next = &(fCounter2->fNext);
2765  while(*next != 0) next = &( (*next)->fNext);
2766  *next = new TFormLeafInfoClones(motherclassptr,offset,element);
2767 }
2768 
2769 ////////////////////////////////////////////////////////////////////////////////
2770 /// Copy the object and all its data.
2771 
2773 {
2774  return new TFormLeafInfoMultiVarDimClones(*this);
2775 }
2776 
2777 ////////////////////////////////////////////////////////////////////////////////
2778 
2780  Int_t /* instance */)
2781 {
2782  /* The proper indexing and unwinding of index need to be done by prior leafinfo in the chain. */
2783  Error("GetValue","This should never be called");
2784  return 0;
2785 }
2786 
2787 ////////////////////////////////////////////////////////////////////////////////
2788 /// Load the current array sizes.
2789 
2791 {
2793 
2794  TLeaf *leaf = (TLeaf*)branch->GetListOfLeaves()->At(0);
2796 
2797  if (fNsize > fSizes.GetSize()) fSizes.Set(fNsize);
2798  fSumOfSizes = 0;
2799  for (Int_t i=0; i<fNsize; i++) {
2800  TClonesArray *clones = (TClonesArray*)fCounter2->GetValuePointer(leaf,i);
2801  if (clones) {
2802  Int_t size = clones->GetEntries();
2803  fSumOfSizes += size;
2804  fSizes.AddAt( size, i );
2805  }
2806  }
2807  return;
2808 }
2809 
2810 ////////////////////////////////////////////////////////////////////////////////
2811 /// Return the value of the underlying data.
2812 
2813 template <typename T>
2814 T TFormLeafInfoMultiVarDimClones::ReadValueImpl(char *where, Int_t instance)
2815 {
2816  if (fSecondaryIndex>=0) {
2817  UInt_t len = fNext->GetArrayLength();
2818  if (len) {
2819  instance = fSecondaryIndex*len;
2820  } else {
2821  instance = fSecondaryIndex;
2822  }
2823  }
2824  return fNext->ReadTypedValue<T>(where,instance);
2825 }
2826 
2828 
2829 ////////////////////////////////////////////////////////////////////////////////
2830 /// \class TFormLeafInfoCast
2831 /// A small helper class to implement casting an object to
2832 /// a different type (equivalent to dynamic_cast)
2833 
2834 ////////////////////////////////////////////////////////////////////////////////
2835 /// Constructor.
2836 
2838  TFormLeafInfo(classptr),fCasted(casted),fGoodCast(kTRUE)
2839 {
2840  if (casted) { fCastedName = casted->GetName(); }
2841  fMultiplicity = -1;
2843 }
2844 
2845 ////////////////////////////////////////////////////////////////////////////////
2846 /// Constructor.
2847 
2849  TFormLeafInfo(orig)
2850 {
2851  fCasted = orig.fCasted;
2852  fCastedName = orig.fCastedName;
2853  fGoodCast = orig.fGoodCast;
2854  fIsTObject = orig.fIsTObject;
2855 }
2856 
2857 ////////////////////////////////////////////////////////////////////////////////
2858 /// Exception safe swap.
2859 
2861 {
2862  TFormLeafInfo::Swap(other);
2863  std::swap(fCasted,other.fCasted);
2865  std::swap(fGoodCast,other.fGoodCast);
2867 }
2868 
2869 ////////////////////////////////////////////////////////////////////////////////
2870 /// Exception safe assignment operator.
2871 
2873 {
2874  TFormLeafInfoCast tmp(other);
2875  Swap(tmp);
2876  return *this;
2877 }
2878 
2879 ////////////////////////////////////////////////////////////////////////////////
2880 /// Copy the object and all its content.
2881 
2883 {
2884  return new TFormLeafInfoCast(*this);
2885 }
2886 
2887 ////////////////////////////////////////////////////////////////////////////////
2888 /// Destructor.
2889 
2891 {
2892 }
2893 
2894 // Currently only implemented in TFormLeafInfoCast
2896 {
2897  // Get the number of element in the entry.
2898 
2899  if (!fGoodCast) return 0;
2900  if (fNext) return fNext->GetNdata();
2901  return 1;
2902 }
2903 
2904 ////////////////////////////////////////////////////////////////////////////////
2905 /// Read the value at the given memory location
2906 
2907 template <typename T>
2908 T TFormLeafInfoCast::ReadValueImpl(char *where, Int_t instance)
2909 {
2910  if (!fNext) return 0;
2911 
2912  // First check that the real class inherits from the
2913  // casted class
2914  // First assume TObject ...
2915  if ( fIsTObject && !((TObject*)where)->InheritsFrom(fCasted) ) {
2916  fGoodCast = kFALSE;
2917  return 0;
2918  } else {
2919  // We know we have a TBranchElement and we need to find out the
2920  // real class name.
2921  }
2922  fGoodCast = kTRUE;
2923  return fNext->ReadTypedValue<T>(where,instance);
2924 }
2925 
2927 
2928 
2929 ////////////////////////////////////////////////////////////////////////////////
2930 /// We reloading all cached information in case the underlying class
2931 /// information has changed (for example when changing from the 'emulated'
2932 /// class to the real class.
2933 
2935 {
2936  if (fCasted) {
2937  TClass * new_class = TClass::GetClass(fCastedName);
2938  if (new_class!=fCasted) {
2939  fCasted = new_class;
2940  }
2941  }
2942  return TFormLeafInfo::Update();
2943 }
2944 
2945 ////////////////////////////////////////////////////////////////////////////////
2946 /// \class TFormLeafInfoTTree
2947 /// A small helper class to implement reading
2948 /// from the containing TTree object itself.
2949 
2950 TFormLeafInfoTTree::TFormLeafInfoTTree(TTree *tree, const char *alias, TTree *current) :
2951 TFormLeafInfo( TTree::Class(), 0, 0 ), fTree(tree),fCurrent(current),fAlias(alias)
2952 {
2953  // Constructor.
2954 
2955  if (fCurrent==0) fCurrent = fTree->GetFriend(alias);
2956 }
2957 
2959  TFormLeafInfo(orig)
2960 {
2961  // Copy Constructor.
2962  fTree = orig.fTree;
2963  fAlias = orig.fAlias;
2964  fCurrent = orig.fCurrent;
2965 }
2966 
2968 {
2969  // Copy the object and all its content.
2970 
2971  return new TFormLeafInfoTTree(*this);
2972 }
2973 
2974 ////////////////////////////////////////////////////////////////////////////////
2975 /// returns the address of the value pointed to by the
2976 /// TFormLeafInfo.
2977 
2979 {
2980  return GetLocalValuePointer((char*)fCurrent,instance);
2981 }
2982 
2983 ////////////////////////////////////////////////////////////////////////////////
2984 /// Return result of a leafobject method.
2985 
2986 template <typename T>
2987 T TFormLeafInfoTTree::GetValueImpl(TLeaf *, Int_t instance)
2988 {
2989  return ReadTypedValue<T>((char*)fCurrent,instance);
2990 }
2991 
2992 ////////////////////////////////////////////////////////////////////////////////
2993 /// Return result of a leafobject method.
2994 
2995 template <typename T>
2996 T TFormLeafInfoTTree::ReadValueImpl(char *thisobj, Int_t instance)
2997 {
2998  if (fElement) return TFormLeafInfo::ReadTypedValue<T>(thisobj,instance);
2999  else if (fNext) return fNext->ReadTypedValue<T>(thisobj,instance);
3000  else return 0;
3001 }
3002 
3005 
3006 ////////////////////////////////////////////////////////////////////////////////
3007 /// Update after a change of file in a chain
3008 
3010 {
3011  if (fAlias.Length() && fAlias != fTree->GetName()) {
3013  }
3014  return fCurrent && TFormLeafInfo::Update();
3015 }
Describe Streamer information for one class version.
Definition: TStreamerInfo.h:47
Int_t GetNdata() const
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:37
virtual Int_t ReadCounterValue(char *where)
Return the current size of the the TClonesArray.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
TFormLeafInfoNumerical & operator=(const TFormLeafInfoNumerical &orig)
Exception safe assignment operator.
#define INSTANTIATE_GETVAL(CLASS)
TFormLeafInfo & operator=(const TFormLeafInfo &orig)
Exception safe assignment operator.
~TFormLeafInfoCollection()
Destructor.
A small helper class to implement reading a data member on a TClonesArray object stored in a TTree...
virtual void AddOffset(Int_t offset, TStreamerElement *element)
Increase the offset of this element.
long long Long64_t
Definition: RtypesCore.h:69
Int_t GetID() const
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
TVirtualCollectionProxy * fCollProxy
Bool_t fGoodCast
Name of the class we are casting to.
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:488
Int_t GetNdata(TLeaf *leaf)
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
Ssiz_t Length() const
Definition: TString.h:390
TFormLeafInfoCollectionSize()
Constructor.
float Float_t
Definition: RtypesCore.h:53
Equal to TDataType&#39;s kchar.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
static const EReturnType kOther
Definition: TMethodCall.h:50
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
void swap(TDirectoryEntry &e1, TDirectoryEntry &e2) noexcept
A small helper class to implement reading a data member on an object stored in a TTree.
double T(double x)
Definition: ChebyshevPol.h:34
TString fClassName
Definition: TFormLeafInfo.h:74
TBranchElement * GetBranchCount2() const
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
returns the address of the value pointed to by the serie of TFormLeafInfo.
A TLeaf for a general object derived from TObject.
Definition: TLeafObject.h:35
char * GetObject() const
Return a pointer to our object.
const char * GetParams() const
Definition: TMethodCall.h:95
unsigned short UShort_t
Definition: RtypesCore.h:36
const char * GetReturnTypeName() const
Get full type description of function return type, e,g.: "class TDirectory*".
Definition: TFunction.cxx:140
Bool_t IsTObject() const
Return kTRUE is the class inherits from TObject.
Definition: TClass.cxx:5571
virtual Double_t GetValue(TLeaf *leaf, Int_t instance=0)
TFormLeafInfoCollectionObject(TClass *classptr=0, Bool_t fTop=kTRUE)
Constructor.
void Swap(TFormLeafInfoCast &other)
Exception safe swap.
virtual EDataType GetType() const =0
#define R__ASSERT(e)
Definition: TError.h:98
#define gROOT
Definition: TROOT.h:364
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all of its content.
void Swap(TFormLeafInfoMultiVarDim &other)
Exception safe swap.
static const EReturnType kLong
Definition: TMethodCall.h:47
A small helper class to implement casting an object to a different type (equivalent to dynamic_cast) ...
Basic string class.
Definition: TString.h:137
virtual TTree * GetFriend(const char *) const
Return a pointer to the TTree friend whose name or alias is &#39;friendname.
Definition: TTree.cxx:5496
TFormLeafInfoCast & operator=(const TFormLeafInfoCast &orig)
Exception safe assignment operator.
TFormLeafInfo * fNext
Definition: TFormLeafInfo.h:73
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
#define gInterpreter
Definition: TInterpreter.h:515
TBranch * GetBranch() const
Definition: TLeaf.h:70
TFormLeafInfoNumerical(TVirtualCollectionProxy *holder_of)
Construct a TFormLeafInfo for the numerical type contained in the collection.
TFormLeafInfo(TClass *classptr=0, Long_t offset=0, TStreamerElement *element=0)
Constructor.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
virtual Int_t GetSize(Int_t index)
Return the size of the requested sub-array.
const char * GetMethodName() const
Definition: TMethodCall.h:94
const char * Class
Definition: TXMLSetup.cxx:64
~TFormLeafInfoCollectionSize()
Destructor.
virtual Int_t GetPrimaryIndex()
Return the index vlaue of the primary index.
TFormLeafInfoClones(TClass *classptr=0, Long_t offset=0)
Constructor.
virtual void * GetLocalValuePointer(TLeaf *from, Int_t instance=0)
This is implemented here because some compiler want ALL the signature of an overloaded function to be...
virtual TFormLeafInfo * DeepCopy() const
Make a complete copy of this FormLeafInfo and all its content.
TFormLeafInfoCollection()
Constructor.
TFormLeafInfoCollectionObject & operator=(const TFormLeafInfoCollectionObject &orig)
Exception safe assignement operator.
#define READ_ARRAY(TYPE_t)
A helper class to implement reading a data member on a variable size array inside a TClonesArray obje...
Array of integers (32 bits per element).
Definition: TArrayI.h:29
virtual TClass * GetCollectionClass() const
TMethodCall * fMethod
static const EReturnType kString
Definition: TMethodCall.h:49
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual TClass * GetClass() const
Get the class of the underlying data.
TBranchElement * GetBranchCount() const
TFormLeafInfoMultiVarDim()
Constructor.
const char * Data() const
Definition: TString.h:349
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all of its contnet.
TFormLeafInfoMethod & operator=(const TFormLeafInfoMethod &orig)
Exception safe assignment operator.
virtual Bool_t IsString() const
Return true if the underlying data is a string.
virtual Int_t GetVarDim()
Return the index of the dimension which varies for each elements of an enclosing array (typically a T...
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual Bool_t HasPointers() const =0
TFormLeafInfoMultiVarDim & operator=(const TFormLeafInfoMultiVarDim &orig)
Exception safe assignment operator.
virtual TFormLeafInfo * DeepCopy() const
Copy this object and its content.
TStreamerInfo * GetInfo() const
Get streamer info for the branch class.
virtual Bool_t Update()
Update after a change of file in a chain.
TString fCastedName
Pointer to the class we are trying to case to.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual Int_t GetSumOfSizes()
Get the total size.
void Init(const TFunction *func)
Initialize the method invocation environment based on the TFunction object.
TString fElementName
Definition: TFormLeafInfo.h:75
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Not implemented.
virtual TVirtualCollectionProxy * Generate() const =0
virtual ~TFormLeafInfoCast()
Destructor.
TFormLeafInfoTTree(TTree *tree=0, const char *alias=0, TTree *current=0)
TFormLeafInfoMultiVarDimCollection(TClass *motherclassptr, Long_t offset, TClass *elementclassptr, TFormLeafInfo *parent)
Constructor.
virtual Double_t ReadValue(char *, Int_t=0)
Read the value at the given memory location.
TString & Append(const char *cs)
Definition: TString.h:492
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
virtual void SetPrimaryIndex(Int_t index)
Set the current value of the primary index.
static const EReturnType kDouble
Definition: TMethodCall.h:48
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
virtual Double_t ReadValue(char *where, Int_t instance=0)
Return the value of the underlying data member inside the clones array.
virtual Bool_t IsInteger() const
Return true if the underlying data is an integral value.
Method or function calling interface.
Definition: TMethodCall.h:41
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
void Set(Int_t n)
Set size of this array to n ints.
Definition: TArrayI.cxx:105
void AddAt(Int_t c, Int_t i)
Add Int_t c at position i. Check for out of bounds.
Definition: TArrayI.cxx:93
virtual void UpdateSizes(TArrayI *garr)
Update the sizes of the arrays.
virtual Int_t ReadCounterValue(char *where)
Return the size of the underlying array for the current entry in the TTree.
TFormLeafInfoMethod(TClass *classptr=0, TMethodCall *method=0)
Constructor.
TObject * UncheckedAt(Int_t i) const
Definition: TObjArray.h:91
char * GetObjectAddress(TLeafElement *leaf, Int_t &instance)
Returns the the location of the object pointed to.
Int_t fMultiplicity
Definition: TFormLeafInfo.h:78
A small helper class to implement reading a data member on a variable size array inside a TClonesArra...
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the current size of the the TClonesArray.
Used to return the size of a collection.
TClass * fClass
Definition: TFormLeafInfo.h:66
virtual ~TFormLeafInfo()
Delete this object and all its content.
virtual void SetSecondaryIndex(Int_t index)
Set the primary index value.
TRandom2 r(17)
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
A small helper class to implement reading a data member on a variable size array inside a TClonesArra...
TVirtualStreamerInfo * GetStreamerInfo(Int_t version=0) const
returns a pointer to the TVirtualStreamerInfo object for version If the object does not exist...
Definition: TClass.cxx:4356
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the current size of the the TClonesArray.
T GetTypedValue(Int_t i=0) const
Definition: TLeaf.h:86
TStreamerElement * fElement
Offset of the data pointed inside the class fClass.
Definition: TFormLeafInfo.h:69
A small helper class to implement reading a data member by following a pointer inside a branch of TTr...
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the current size of the the TClonesArray.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual void SetSize(Int_t index, Int_t val)
Set the current size of the arrays.
void Copy(TArrayI &array) const
Definition: TArrayI.h:44
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
virtual Int_t GetPrimaryIndex()
Method for multiple variable dimensions.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
returns the address of the value pointed to by the TFormLeafInfo.
virtual void SetSize(Int_t index, Int_t val)
Set the sizes of the sub-array.
virtual Bool_t HasCounter() const
Return true if the underlying data has a array size counter.
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
short Short_t
Definition: RtypesCore.h:35
TLine * l
Definition: textangle.C:4
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:81
T ReadTypedValue(char *where, Int_t instance=0)
virtual Int_t GetSize(Int_t index)
For the current entry, and the value &#39;index&#39; for the main array, return the size of the secondary var...
Long_t Property() const
Get property description word. For meaning of bits see EProperty.
Definition: TFunction.cxx:183
virtual void UpdateSizes(TArrayI *garr)
Set the current sizes of the arrays.
virtual void SetPrimaryIndex(Int_t index)
Set the primary index value.
TVirtualCollectionProxy * fCollProxy
#define INSTANTIATE_READVAL(CLASS)
TFormLeafInfoClones & operator=(const TFormLeafInfoClones &orig)
Exception safe assignement operator.
virtual ~TFormLeafInfoNumerical()
Destructor.
virtual Double_t ReadValue(char *, Int_t=0)
Not implemented.
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Not implemented.
TTree * GetTree() const
Definition: TBranch.h:184
A TLeaf for the general case when using the branches created via a TStreamerInfo (i.e.
Definition: TLeafElement.h:34
A Branch for the case of an object.
virtual void SetSecondaryIndex(Int_t index)
Set the current value of the primary index.
Int_t GetSize() const
Definition: TArray.h:49
TFormLeafInfoCollection & operator=(const TFormLeafInfoCollection &orig)
Exception safe assignment operator.
TFormLeafInfoCollectionSize & operator=(const TFormLeafInfoCollectionSize &orig)
Exception safe assignment operator.
~TFormLeafInfoMethod()
Destructor.
void Swap(TFormLeafInfoNumerical &other)
Exception safe swap.
long Long_t
Definition: RtypesCore.h:50
Bool_t IsLoaded() const
Return true if the shared library of this class is currently in the a process&#39;s memory.
Definition: TClass.cxx:5545
virtual TFormLeafInfo * DeepCopy() const
Make a complete copy of this FormLeafInfo and all its content.
virtual Int_t GetVirtVarDim()
Return the virtual index (for this expression) of the dimension which varies for each elements of an ...
TFormLeafInfoPointer(TClass *classptr=0, Long_t offset=0, TStreamerElement *element=0)
Constructor.
const char * GetCountName() const
TFunction * GetMethod()
Returns the TMethod describing the method to be executed.
A small helper class to implement reading from the containing TTree object itself.
virtual TFormLeafInfo * DeepCopy() const
Copy of the object and its content.
virtual Int_t ReadCounterValue(char *where)
Return the size of the underlying array for the current entry in the TTree.
virtual Int_t GetNdata()
Get the number of element in the entry.
TFormLeafInfoMultiVarDimClones(TClass *motherclassptr, Long_t offset, TClass *elementclassptr, TFormLeafInfo *parent)
Constructor.
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:416
double Double_t
Definition: RtypesCore.h:55
virtual TClass * GetClassPointer() const
Returns a pointer to the TClass of this element.
A small helper class to implement reading a numerical value inside a collection.
virtual Bool_t IsString() const
Return true if the return value is a string.
Int_t At(Int_t i) const
Definition: TArrayI.h:81
This class is a small helper class to implement reading a data member on an object stored in a TTree...
Definition: TFormLeafInfo.h:53
int type
Definition: TGX11.cxx:120
Int_t GetNewType(Int_t id) const
unsigned long long ULong64_t
Definition: RtypesCore.h:70
Int_t GetMultiplicity()
Reminder of the meaning of fMultiplicity:
unsigned long ULong_t
Definition: RtypesCore.h:51
EDataType
Definition: TDataType.h:30
TFormLeafInfoDirect(TBranchElement *from)
Constructor.
Int_t GetNewType() const
virtual void * At(UInt_t idx)=0
virtual Int_t GetNdata()
Get the number of element in the entry.
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
TObjArray * GetListOfLeaves()
Definition: TBranch.h:179
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:494
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition: TClass.cxx:2882
virtual Int_t GetSumOfSizes()
Total all the elements that are available for the current entry for the secondary variable dimension...
void Swap(TFormLeafInfoCollectionSize &other)
Exception safe swap.
TFormLeafInfo * fCounter2
TVirtualCollectionProxy * GetCollectionProxy() const
Return the proxy describing the collection (if any).
Definition: TClass.cxx:2811
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
void Swap(TFormLeafInfoCollection &other)
Exception safe swap.
TClass * GetClass() const
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
Mother of all ROOT objects.
Definition: TObject.h:44
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the address of the underlying value.
char Char_t
Definition: RtypesCore.h:29
void Swap(TFormLeafInfoClones &other)
Exception safe swap.
virtual Double_t ReadValue(char *where, Int_t instance=0)
Return the value of the underlying pointer data member.
void Swap(TFormLeafInfo &other)
virtual Double_t GetValue(TLeaf *, Int_t=0)
An array of clone (identical) objects.
Definition: TClonesArray.h:32
virtual Bool_t HasCounter() const
Return true if any of underlying data has a array size counter.
virtual char * GetAddress() const
Get the branch address.
virtual void * GetValuePointer() const
Definition: TLeaf.h:80
void Execute(const char *, const char *, int *=0)
Execute method on this object with the given parameter string, e.g.
Definition: TMethodCall.h:68
TFormLeafInfoCast(TClass *classptr=0, TClass *casted=0)
Indicated whether the fClass inherits from TObject.
virtual Bool_t IsInteger() const
Return true if the return value is integral.
TStreamerElement * fLocalElement
Definition: tree.py:1
A small helper class to implement reading a data member on a generic collection object stored in a TT...
A TTree object has a header with a name and a title.
Definition: TTree.h:98
Asmall helper class to implement executing a method of an object stored in a TTree.
double result[121]
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
returns the address of the value pointed to by the TFormLeafInfo.
unsigned char UChar_t
Definition: RtypesCore.h:34
Int_t GetMakeClass() const
Definition: TTree.h:423
Bool_t fIsTObject
Marked by ReadValue.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
TFormLeafInfo * fCounter
Descriptor of the data pointed to.
Definition: TFormLeafInfo.h:72
Long_t fOffset
This is the class of the data pointed to.
Definition: TFormLeafInfo.h:68
TObject * At(Int_t idx) const
Definition: TObjArray.h:167
virtual Int_t GetArrayLength()
Return the current length of the array.
virtual Bool_t IsString() const
Return true if the underlying data is a string.
A TTree is a list of TBranches.
Definition: TBranch.h:58
A small helper class to implement reading a data member on a TClonesArray object stored in a TTree...
const Bool_t kTRUE
Definition: Rtypes.h:91
void Swap(TFormLeafInfoCollectionObject &other)
Exception safe swap.
A small helper class to implement reading a data member on a variable size array inside a TClonesArra...
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the size of the underlying array for the current entry in the TTree.
virtual TClass * GetClass() const
Return the type of the underlying return value.
virtual Int_t GetVarDim()
Return the index of the dimension which varies for each elements of an enclosing array (typically a T...
virtual Int_t GetVirtVarDim()
Return the virtual index (for this expression) of the dimension which varies for each elements of an ...
Int_t Size() const
Return size of object of this class.
Definition: TClass.cxx:5344
EReturnType ReturnType()
Returns the return type of the method.
~TFormLeafInfoMultiVarDim()
Destructor.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its data.
#define GET_ARRAY(TYPE_t)
virtual UInt_t Size() const =0
Int_t GetArrayLength() const
void Swap(TFormLeafInfoMethod &other)
Exception safe swap.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
virtual Double_t GetValue(TLeaf *leaf, Int_t instance=0)