Logo ROOT  
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
13This class is a small helper class to implement reading a data member
14on an object stored in a TTree.
15
16TTreeFormula now relies on a variety of TFormLeafInfo classes to handle the
17reading 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
29The 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#include "TClassEdit.h"
58
59
60#define INSTANTIATE_READVAL(CLASS) \
61 template Double_t CLASS::ReadValueImpl<Double_t>(char*, Int_t); \
62 template Long64_t CLASS::ReadValueImpl<Long64_t>(char*, Int_t); \
63 template LongDouble_t CLASS::ReadValueImpl<LongDouble_t>(char*, Int_t) // no semicolon
64
65
66#define INSTANTIATE_GETVAL(CLASS) \
67 template Double_t CLASS::GetValueImpl<Double_t>(TLeaf*, Int_t); \
68 template Long64_t CLASS::GetValueImpl<Long64_t>(TLeaf*, Int_t); \
69 template LongDouble_t CLASS::GetValueImpl<LongDouble_t>(TLeaf*, Int_t) // no semicolon
70
71////////////////////////////////////////////////////////////////////////////////
72/// Constructor.
73
75 TStreamerElement* element) :
76 fClass(classptr),fOffset(offset),fElement(element),
77 fCounter(0), fNext(0),fMultiplicity(0)
78{
80 if (fElement) {
82 }
83}
84
85////////////////////////////////////////////////////////////////////////////////
86///Constructor.
87
88TFormLeafInfo::TFormLeafInfo(const TFormLeafInfo& orig) : TObject(orig),fClass(orig.fClass),fOffset(orig.fOffset),fElement(orig.fElement),fCounter(0),fNext(0),fClassName(orig.fClassName),fElementName(orig.fElementName),fMultiplicity(orig.fMultiplicity)
89{
90 // Deep copy the pointers.
91 if (orig.fCounter) fCounter = orig.fCounter->DeepCopy();
92 if (orig.fNext) fNext = orig.fNext->DeepCopy();
93}
95////////////////////////////////////////////////////////////////////////////////
96/// Exception safe assignment operator.
97
99{
100 TFormLeafInfo tmp(other);
101 Swap(tmp);
102 return *this;
103}
104
105////////////////////////////////////////////////////////////////////////////////
106
108{
109 std::swap(fClass,other.fClass);
113 std::swap(fNext,other.fNext);
114 TString tmp(fClassName);
115 fClassName = other.fClassName;
116 other.fClassName = tmp;
117
118 tmp = fElementName;
120 other.fElementName = tmp;
122}
123
124////////////////////////////////////////////////////////////////////////////////
125/// Make a complete copy of this FormLeafInfo and all its content.
126
128{
129 return new TFormLeafInfo(*this);
130}
131
132////////////////////////////////////////////////////////////////////////////////
133/// Delete this object and all its content
134
136{
137 delete fCounter;
138 delete fNext;
139}
140
141////////////////////////////////////////////////////////////////////////////////
142/// Increase the offset of this element. This intended to be the offset
143/// from the start of the object to which the data member belongs.
144
146{
147 fOffset += offset;
148 fElement = element;
149 if (fElement ) {
150 // fElementClassOwnerName = cl->GetName();
151 fElementName.Append(".").Append(element->GetName());
152 }
153}
154
155////////////////////////////////////////////////////////////////////////////////
156/// Return the current length of the array.
157
159{
160 Int_t len = 1;
161 if (fNext) len = fNext->GetArrayLength();
162 if (fElement) {
163 Int_t elen = fElement->GetArrayLength();
164 if (elen || fElement->IsA() == TStreamerBasicPointer::Class() )
165 len *= fElement->GetArrayLength();
166 }
167 return len;
168}
169
170////////////////////////////////////////////////////////////////////////////////
171/// Get the class of the underlying data.
172
174{
175 if (fNext) return fNext->GetClass();
176 if (fElement) return fElement->GetClassPointer();
177 return fClass;
178}
179
180////////////////////////////////////////////////////////////////////////////////
181/// Returns the the location of the object pointed to.
182/// Modify instance if the object is part of an array.
183
185{
186 TBranchElement* branch = (TBranchElement*) leaf->GetBranch();
187 Int_t id = branch->GetID();
188 if (id < 0) {
189 // Branch is a top-level branch.
190 if (branch->GetTree()->GetMakeClass()) {
191 // Branch belongs to a MakeClass tree.
192 return branch->GetAddress();
193 } else {
194 return branch->GetObject();
195 }
196 }
197 TStreamerInfo* info = branch->GetInfo();
198 Int_t offset = 0;
199 if (id > -1) {
200 // Branch is *not* a top-level branch.
201 offset = info->TStreamerInfo::GetElementOffset(id);
202 }
203 char* address = 0;
204 // Branch is *not* a top-level branch.
205 if (branch->GetTree()->GetMakeClass()) {
206 // Branch belongs to a MakeClass tree.
207 address = (char*) branch->GetAddress();
208 } else {
209 address = (char*) branch->GetObject();
210 }
211 char* thisobj = 0;
212 if (!address) {
213 // FIXME: This makes no sense, if the branch address is not set, then object will not be set either.
214 thisobj = branch->GetObject();
215 } else {
216 Int_t type = -1;
217 if (id > -1) {
218 // Note this is somewhat slow
219 type = info->TStreamerInfo::GetElement(id)->GetNewType();
220 }
221 switch (type) {
227 Error("GetValuePointer", "Type (%d) not yet supported\n", type);
228 break;
229
233 {
234 // An array of objects.
235 Int_t index;
236 Int_t sub_instance;
237 Int_t len = GetArrayLength();
238 if (len) {
239 index = instance / len;
240 sub_instance = instance % len;
241 } else {
242 index = instance;
243 sub_instance = 0;
244 }
245 thisobj = address + offset + (index * fClass->Size());
246 instance = sub_instance;
247 break;
248 }
249
257 // A single object.
258 thisobj = address + offset;
259 break;
260
294 // A simple type, or an array of a simple type.
295 thisobj = address + offset;
296 break;
297
298 default:
299 // Everything else is a pointer to something.
300 thisobj = *((char**) (address + offset));
301 break;
302 }
303 }
304 return thisobj;
305}
306
307////////////////////////////////////////////////////////////////////////////////
308/// Reminder of the meaning of fMultiplicity:
309/// - -1: Only one or 0 element per entry but contains variable length array!
310/// - 0: Only one element per entry, no variable length array
311/// - 1: loop over the elements of a variable length array
312/// - 2: loop over elements of fixed length array (nData is the same for all entry)
313
315{
316 // Currently only TFormLeafInfoCast uses this field.
317 return fMultiplicity;
318}
319
321{
322 // Get the number of element in the entry.
323
324 GetCounterValue(leaf);
325 GetValue(leaf);
326 return GetNdata();
327}
328
329////////////////////////////////////////////////////////////////////////////////
330/// Get the number of element in the entry.
331
333{
334 if (fNext) return fNext->GetNdata();
335 return 1;
336}
337
338////////////////////////////////////////////////////////////////////////////////
339/// Return true if any of underlying data has a array size counter
340
342{
343 Bool_t result = kFALSE;
344 if (fNext) result = fNext->HasCounter();
345 return fCounter!=0 || result;
346}
347
348////////////////////////////////////////////////////////////////////////////////
349/// Return true if the underlying data is a string
350
352{
353 if (fNext) return fNext->IsString();
354 if (!fElement) return kFALSE;
355
356 switch (fElement->GetNewType()) {
357 // basic types
358 case kChar_t:
359 // This is new in ROOT 3.02/05
360 return kFALSE;
362 // This is new in ROOT 3.02/05
363 return kTRUE;
365 return kTRUE;
366 default:
367 return kFALSE;
368 }
369}
370
371////////////////////////////////////////////////////////////////////////////////
372/// Return true if the underlying data is an integral value
373
375{
376 if (fNext) return fNext->IsInteger();
377 if (!fElement) return kFALSE;
378
379 Int_t atype = fElement->GetNewType();
380 if (TStreamerInfo::kOffsetL < atype &&
381 atype < TStreamerInfo::kOffsetP ) {
383 } else if (TStreamerInfo::kOffsetP < atype &&
384 atype < TStreamerInfo::kObject) {
386 }
387
388 switch (atype) {
389 // basic types
402 return kTRUE;
404 return kTRUE; // For consistency with the leaf list method and proper axis setting
409 return kFALSE;
410 default:
411 return kFALSE;
412 }
413}
414
415////////////////////////////////////////////////////////////////////////////////
416/// Method for multiple variable dimensions.
417
419{
420 if (fNext) return fNext->GetPrimaryIndex();
421 return -1;
422}
423
424////////////////////////////////////////////////////////////////////////////////
425/// Return the index of the dimension which varies
426/// for each elements of an enclosing array (typically a TClonesArray)
427
429{
430 if (fNext) return fNext->GetVarDim();
431 else return -1;
432}
433
434////////////////////////////////////////////////////////////////////////////////
435/// Return the virtual index (for this expression) of the dimension which varies
436/// for each elements of an enclosing array (typically a TClonesArray)
437
439{
440 if (fNext) return fNext->GetVirtVarDim();
441 else return -1;
442}
443
444////////////////////////////////////////////////////////////////////////////////
445/// For the current entry, and the value 'index' for the main array,
446/// return the size of the secondary variable dimension of the 'array'.
447
449{
450 if (fNext) return fNext->GetSize(index);
451 else return 0;
452}
453
454////////////////////////////////////////////////////////////////////////////////
455/// Total all the elements that are available for the current entry
456/// for the secondary variable dimension.
457
459{
460 if (fNext) return fNext->GetSumOfSizes();
461 else return 0;
462}
463
464////////////////////////////////////////////////////////////////////////////////
465/// Load the current array sizes
466
468{
469 if (fNext) fNext->LoadSizes(branch);
470}
471
472////////////////////////////////////////////////////////////////////////////////
473/// Set the primary index value
474
476{
477 if (fNext) fNext->SetPrimaryIndex(index);
478}
479
480////////////////////////////////////////////////////////////////////////////////
481/// Set the primary index value
482
484{
485 if (fNext) fNext->SetSecondaryIndex(index);
486}
487
488////////////////////////////////////////////////////////////////////////////////
489/// Set the current size of the arrays
490
492{
493 if (fNext) fNext->SetSize(index, val);
494}
495
496////////////////////////////////////////////////////////////////////////////////
497/// Set the current sizes of the arrays
498
500{
501 if (fNext) fNext->UpdateSizes(garr);
502}
503
504////////////////////////////////////////////////////////////////////////////////
505/// We reloading all cached information in case the underlying class
506/// information has changed (for example when changing from the 'emulated'
507/// class to the real class.
508
510{
511 if (fClass) {
512 TClass * new_class = TClass::GetClass(fClassName);
513 if (new_class==fClass) {
514 if (fNext) fNext->Update();
515 if (fCounter) fCounter->Update();
516 return kFALSE;
517 }
518 fClass = new_class;
519 }
520 if (fElement && fClass) {
521 TClass *cl = fClass;
522 // We have to drill down the element name within the class.
523 Int_t offset,i;
524 TStreamerElement* element;
525 char * current;
526 Int_t nchname = fElementName.Length();
527 char * work = new char[nchname+2];
528 for (i=0, current = &(work[0]), fOffset=0; i<nchname+1;i++ ) {
529 if (i==nchname || fElementName[i]=='.') {
530 // A delimiter happened let's see if what we have seen
531 // so far does point to a data member.
532 *current = '\0';
533 element = ((TStreamerInfo*)cl->GetStreamerInfo())->GetStreamerElement(work,offset);
534 if (element) {
535 Int_t type = element->GetNewType();
536 if (type<60) {
537 fOffset += offset;
538 } else if (type == TStreamerInfo::kBase ||
555 fOffset += offset;
556 cl = element->GetClassPointer();
557 }
558 fElement = element;
559 current = &(work[0]);
560 }
561 } else {
562 if (i<nchname) *current++ = fElementName[i];
563 }
564 }
565 delete [] work;
566 }
567 if (fNext) fNext->Update();
568 if (fCounter) fCounter->Update();
569 return kTRUE;
570}
571
572////////////////////////////////////////////////////////////////////////////////
573/// Return the size of the underlying array for the current entry in the TTree.
574
576 if (!fCounter) {
577 if (fNext && fNext->HasCounter()) {
578 char *where = (char*)GetLocalValuePointer(leaf,0);
579 return fNext->ReadCounterValue(where);
580 } else return 1;
581 }
582 return (Int_t)fCounter->GetValue(leaf);
583}
584
585////////////////////////////////////////////////////////////////////////////////
586/// Return the size of the underlying array for the current entry in the TTree.
587
589{
590 if (!fCounter) {
591 if (fNext) {
592 char *next = (char*)GetLocalValuePointer(where,0);
593 return fNext->ReadCounterValue(next);
594 } else return 1;
595 }
596 return (Int_t)fCounter->ReadValue(where,0);
597}
598
599////////////////////////////////////////////////////////////////////////////////
600/// returns the address of the value pointed to by the
601/// TFormLeafInfo.
602
604{
605 char *thisobj = 0;
606 if (leaf->InheritsFrom(TLeafObject::Class()) ) {
607 thisobj = (char*)((TLeafObject*)leaf)->GetObject();
608 } else {
609 thisobj = GetObjectAddress((TLeafElement*)leaf, instance); // instance might be modified
610 }
611 if (!thisobj) return 0;
612 return GetLocalValuePointer(thisobj, instance);
613}
614
615////////////////////////////////////////////////////////////////////////////////
616/// returns the address of the value pointed to by the
617/// serie of TFormLeafInfo.
618
620{
621 char *thisobj = (char*)GetLocalValuePointer(leaf,instance);
622 if (fNext) return fNext->GetValuePointer(thisobj,instance);
623 else return thisobj;
624}
625
626////////////////////////////////////////////////////////////////////////////////
627/// returns the address of the value pointed to by the
628/// TFormLeafInfo.
629
631{
632 char *where = (char*)GetLocalValuePointer(thisobj,instance);
633 if (fNext) return fNext->GetValuePointer(where,instance);
634 else return where;
635}
636
637////////////////////////////////////////////////////////////////////////////////
638/// returns the address of the value pointed to by the
639/// TFormLeafInfo.
640
642{
643 if (fElement==0 || thisobj==0) return thisobj;
644
645 switch (fElement->GetNewType()) {
646 // basic types
664 return (Int_t*)(thisobj+fOffset);
665
666 // array of basic types array[8]
668 {Bool_t *val = (Bool_t*)(thisobj+fOffset); return &(val[instance]);}
670 {Char_t *val = (Char_t*)(thisobj+fOffset); return &(val[instance]);}
672 {Short_t *val = (Short_t*)(thisobj+fOffset); return &(val[instance]);}
674 {Int_t *val = (Int_t*)(thisobj+fOffset); return &(val[instance]);}
676 {Long_t *val = (Long_t*)(thisobj+fOffset); return &(val[instance]);}
678 {Long64_t *val = (Long64_t*)(thisobj+fOffset); return &(val[instance]);}
680 {Float_t *val = (Float_t*)(thisobj+fOffset); return &(val[instance]);}
682 {Float_t *val = (Float_t*)(thisobj+fOffset); return &(val[instance]);}
684 {Double_t *val = (Double_t*)(thisobj+fOffset); return &(val[instance]);}
686 {Double_t *val = (Double_t*)(thisobj+fOffset); return &(val[instance]);}
688 {UChar_t *val = (UChar_t*)(thisobj+fOffset); return &(val[instance]);}
690 {UShort_t *val = (UShort_t*)(thisobj+fOffset); return &(val[instance]);}
692 {UInt_t *val = (UInt_t*)(thisobj+fOffset); return &(val[instance]);}
694 {ULong_t *val = (ULong_t*)(thisobj+fOffset); return &(val[instance]);}
696 {ULong64_t *val = (ULong64_t*)(thisobj+fOffset); return &(val[instance]);}
697
698#define GET_ARRAY(TYPE_t) \
699 { \
700 Int_t len, sub_instance, index; \
701 if (fNext) len = fNext->GetArrayLength(); \
702 else len = 1; \
703 if (len) { \
704 index = instance / len; \
705 sub_instance = instance % len; \
706 } else { \
707 index = instance; \
708 sub_instance = 0; \
709 } \
710 TYPE_t **val = (TYPE_t**)(thisobj+fOffset); \
711 return &((val[sub_instance])[index]); \
712 }
713
714 // pointer to an array of basic types array[n]
730
732 {char **stringp = (char**)(thisobj+fOffset); return *stringp;}
733
739 {TObject **obj = (TObject**)(thisobj+fOffset); return *obj; }
740
748 {TObject *obj = (TObject*)(thisobj+fOffset); return obj; }
749
753 char *loc = thisobj+fOffset;
754
755 Int_t len, index;
756 //Int_t sub_instance;
757
758 if (fNext) len = fNext->GetArrayLength();
759 else len = 1;
760 if (len) {
761 index = instance / len;
762 // sub_instance = instance % len;
763 } else {
764 index = instance;
765 // sub_instance = 0;
766 }
767
768 loc += index*fElement->GetClassPointer()->Size();
769
770 TObject *obj = (TObject*)(loc);
771 return obj;
772 }
773
779 {TObject *obj = (TObject*)(thisobj+fOffset); return obj; }
780
781 case kOther_t:
782 default: return 0;
783 }
784}
785
786
787////////////////////////////////////////////////////////////////////////////////
788/// Return result of a leafobject method.
789
790template <typename T>
791T TFormLeafInfo::GetValueImpl(TLeaf *leaf, Int_t instance)
792{
793 char *thisobj = 0;
794 if (leaf->InheritsFrom(TLeafObject::Class()) ) {
795 thisobj = (char*)((TLeafObject*)leaf)->GetObject();
796 } else {
797 thisobj = GetObjectAddress((TLeafElement*)leaf, instance); // instance might be modified
798 }
799 if (thisobj==0) return 0;
800 return ReadTypedValue<T>(thisobj,instance);
801}
802
805
806////////////////////////////////////////////////////////////////////////////////
807/// Read the value at the given memory location
808
809template <typename T>
810T TFormLeafInfo::ReadValueImpl(char *thisobj, Int_t instance)
811{
812 if ( !thisobj ) {
813 Error("ReadValue","Invalid data address: result will be wrong");
814 return 0.0; // Or should throw exception/print error ?
815 }
816 if (fNext) {
817 char *nextobj = thisobj+fOffset;
818 Int_t sub_instance = instance;
823 Int_t index;
824 Int_t len = fNext->GetArrayLength();
825 if (len) {
826 index = instance / len;
827 sub_instance = instance % len;
828 } else {
829 index = instance;
830 sub_instance = 0;
831 }
832 nextobj += index*fElement->GetClassPointer()->Size();
833 }
834 return fNext->ReadTypedValue<T>(nextobj,sub_instance);
835 }
836 // return fInfo->ReadValue(thisobj+fOffset,fElement->GetNewType(),instance,1);
837 switch (fElement->GetNewType()) {
838 // basic types
839 case TStreamerInfo::kBool: return (T)(*(Bool_t*)(thisobj+fOffset));
840 case TStreamerInfo::kChar: return (T)(*(Char_t*)(thisobj+fOffset));
841 case TStreamerInfo::kUChar: return (T)(*(UChar_t*)(thisobj+fOffset));
842 case TStreamerInfo::kShort: return (T)(*(Short_t*)(thisobj+fOffset));
843 case TStreamerInfo::kUShort: return (T)(*(UShort_t*)(thisobj+fOffset));
844 case TStreamerInfo::kInt: return (T)(*(Int_t*)(thisobj+fOffset));
845 case TStreamerInfo::kUInt: return (T)(*(UInt_t*)(thisobj+fOffset));
846 case TStreamerInfo::kLong: return (T)(*(Long_t*)(thisobj+fOffset));
847 case TStreamerInfo::kULong: return (T)(*(ULong_t*)(thisobj+fOffset));
848 case TStreamerInfo::kLong64: return (T)(*(Long64_t*)(thisobj+fOffset));
849 case TStreamerInfo::kULong64: return (T)(*(Long64_t*)(thisobj+fOffset)); //cannot cast to ULong64_t with VC++6
850 case TStreamerInfo::kFloat: return (T)(*(Float_t*)(thisobj+fOffset));
851 case TStreamerInfo::kFloat16: return (T)(*(Float_t*)(thisobj+fOffset));
852 case TStreamerInfo::kDouble: return (T)(*(Double_t*)(thisobj+fOffset));
853 case TStreamerInfo::kDouble32: return (T)(*(Double_t*)(thisobj+fOffset));
854 case TStreamerInfo::kLegacyChar: return (T)(*(char*)(thisobj+fOffset));
855 case TStreamerInfo::kCounter: return (T)(*(Int_t*)(thisobj+fOffset));
856
857 // array of basic types array[8]
859 {Bool_t *val = (Bool_t*)(thisobj+fOffset); return T(val[instance]);}
861 {Char_t *val = (Char_t*)(thisobj+fOffset); return T(val[instance]);}
863 {Short_t *val = (Short_t*)(thisobj+fOffset); return T(val[instance]);}
865 {Int_t *val = (Int_t*)(thisobj+fOffset); return T(val[instance]);}
867 {Long_t *val = (Long_t*)(thisobj+fOffset); return T(val[instance]);}
869 {Long64_t *val = (Long64_t*)(thisobj+fOffset); return T(val[instance]);}
871 {Float_t *val = (Float_t*)(thisobj+fOffset); return T(val[instance]);}
873 {Float_t *val = (Float_t*)(thisobj+fOffset); return T(val[instance]);}
875 {Double_t *val = (Double_t*)(thisobj+fOffset); return T(val[instance]);}
877 {Double_t *val = (Double_t*)(thisobj+fOffset); return T(val[instance]);}
879 {UChar_t *val = (UChar_t*)(thisobj+fOffset); return T(val[instance]);}
881 {UShort_t *val = (UShort_t*)(thisobj+fOffset); return T(val[instance]);}
883 {UInt_t *val = (UInt_t*)(thisobj+fOffset); return T(val[instance]);}
885 {ULong_t *val = (ULong_t*)(thisobj+fOffset); return T(val[instance]);}
886#if defined(_MSC_VER) && (_MSC_VER <= 1200)
888 {Long64_t *val = (Long64_t*)(thisobj+fOffset); return T(val[instance]);}
889#else
891 {ULong64_t *val = (ULong64_t*)(thisobj+fOffset); return T(val[instance]);}
892#endif
893
894#define READ_ARRAY(TYPE_t) \
895 { \
896 Int_t len, sub_instance, index; \
897 len = GetArrayLength(); \
898 if (len) { \
899 index = instance / len; \
900 sub_instance = instance % len; \
901 } else { \
902 index = instance; \
903 sub_instance = 0; \
904 } \
905 TYPE_t **val =(TYPE_t**)(thisobj+fOffset); \
906 return T((val[sub_instance])[index]); \
907 }
908
909 // pointer to an array of basic types array[n]
924#if defined(_MSC_VER) && (_MSC_VER <= 1200)
926#else
928#endif
929
930 case kOther_t:
931 default: return 0;
932 }
933}
934
935////////////////////////////////////////////////////////////////////////////////
936/// \class TFormLeafInfoDirect
937/// A small helper class to implement reading a data
938/// member on an object stored in a TTree.
939
940////////////////////////////////////////////////////////////////////////////////
941/// Constructor.
942
944 TFormLeafInfo(from->GetInfo()->GetClass(),0,
945 from->GetInfo()->GetElement(from->GetID()))
946{
947}
948
949////////////////////////////////////////////////////////////////////////////////
950/// Copy this object and its content.
951
953{
954 return new TFormLeafInfoDirect(*this);
955}
956
957////////////////////////////////////////////////////////////////////////////////
958/// Read the value at the given memory location
959
960Double_t TFormLeafInfoDirect::ReadValue(char * /*where*/, Int_t /*instance*/)
961{
962 Error("ReadValue","Should not be used in a TFormLeafInfoDirect");
963 return 0;
964}
965
966////////////////////////////////////////////////////////////////////////////////
967/// Return the underlying value.
968
969template <typename T>
970T TFormLeafInfoDirect::GetValueImpl(TLeaf *leaf, Int_t instance)
971{
972 return leaf->GetTypedValue<T>(instance);
973}
974
976
977////////////////////////////////////////////////////////////////////////////////
978/// Return the address of the underlying value.
979
981{
982 if (leaf->IsA() != TLeafElement::Class()) {
983 return leaf->GetValuePointer();
984 } else {
985 return GetObjectAddress((TLeafElement*)leaf, instance); // instance might be modified
986 }
987}
988
989////////////////////////////////////////////////////////////////////////////////
990/// Note this should probably never be executed.
991
993{
995}
996
997////////////////////////////////////////////////////////////////////////////////
998/// \class TFormLeafInfoNumerical
999/// A small helper class to implement reading a numerical value inside a collection
1000
1001////////////////////////////////////////////////////////////////////////////////
1002/// Constructor.
1003
1005 TFormLeafInfo(0,0,0),
1006 fKind(kind), fIsBool(kFALSE)
1007{
1008 fElement = new TStreamerElement("data","in collection", 0, fKind, "");
1009}
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// Construct a TFormLeafInfo for the numerical type contained in the collection.
1013
1015 TFormLeafInfo(0,0,0),
1016 fKind(kNoType_t), fIsBool(kFALSE)
1017{
1018 if (collection) {
1019 fKind = (EDataType)collection->GetType();
1021 // Could be a bool
1022 if (strcmp( collection->GetCollectionClass()->GetName(), "vector<bool>") == 0
1023 || strncmp( collection->GetCollectionClass()->GetName(), "bitset<", strlen("bitset<") ) ==0 ) {
1024 fIsBool = kTRUE;
1025 fKind = (EDataType)18;
1026 }
1027 }
1028 }
1029 fElement = new TStreamerElement("data","in collection", 0, fKind, "");
1030}
1031
1032////////////////////////////////////////////////////////////////////////////////
1033/// Constructor.
1034
1036 TFormLeafInfo(orig),
1037 fKind(orig.fKind), fIsBool(kFALSE)
1038{
1039 fElement = new TStreamerElement("data","in collection", 0, fKind, "");
1040}
1041
1042////////////////////////////////////////////////////////////////////////////////
1043/// Exception safe swap.
1044
1046{
1047 TFormLeafInfo::Swap(other);
1048 std::swap(fKind,other.fKind);
1049 std::swap(fIsBool,other.fIsBool);
1050}
1051
1052////////////////////////////////////////////////////////////////////////////////
1053/// Exception safe assignment operator.
1054
1056{
1057 TFormLeafInfoNumerical tmp(other);
1058 Swap(tmp);
1059 return *this;
1060}
1061
1062////////////////////////////////////////////////////////////////////////////////
1063/// Copy the object and all its content.
1064
1066{
1067 return new TFormLeafInfoNumerical(*this);
1068}
1069
1070////////////////////////////////////////////////////////////////////////////////
1071/// Destructor
1072
1074{
1075 delete fElement;
1076}
1077
1078////////////////////////////////////////////////////////////////////////////////
1079/// Return true if the underlying data is a string
1080
1082{
1083 if (fIsBool) return kFALSE;
1084 return TFormLeafInfo::IsString();
1085}
1086
1087////////////////////////////////////////////////////////////////////////////////
1088/// We reloading all cached information in case the underlying class
1089/// information has changed (for example when changing from the 'emulated'
1090/// class to the real class.
1091
1093{
1094 //R__ASSERT(fNext==0);
1095
1096 if (fCounter) return fCounter->Update();
1097 return kFALSE;
1098}
1099
1100namespace {
1101 TStreamerElement *R__GetFakeClonesElem() {
1102 static TStreamerElement gFakeClonesElem("begin","fake",0,
1104 "TClonesArray");
1105 return &gFakeClonesElem;
1106 }
1107}
1108
1109////////////////////////////////////////////////////////////////////////////////
1110/// \class TFormLeafInfoClones
1111/// A small helper class to implement reading a data member
1112/// on a TClonesArray object stored in a TTree.
1113
1114////////////////////////////////////////////////////////////////////////////////
1115/// Constructor.
1116
1118 TFormLeafInfo(classptr,offset,R__GetFakeClonesElem()),fTop(kFALSE)
1119{
1120}
1121
1122////////////////////////////////////////////////////////////////////////////////
1123/// Constructor.
1124
1126 Bool_t top) :
1127 TFormLeafInfo(classptr,offset,R__GetFakeClonesElem()),fTop(top)
1128{
1129}
1130
1131////////////////////////////////////////////////////////////////////////////////
1132/// Constructor.
1133
1135 TStreamerElement* element,
1136 Bool_t top) :
1137 TFormLeafInfo(classptr,offset,element),fTop(top)
1138{
1139}
1140
1141////////////////////////////////////////////////////////////////////////////////
1142/// Deep Copy constructor.
1143
1145 TFormLeafInfo(orig), fTop(orig.fTop)
1146{
1147}
1148
1149////////////////////////////////////////////////////////////////////////////////
1150/// Exception safe swap.
1151
1153{
1154 TFormLeafInfo::Swap(other);
1155 std::swap(fTop,other.fTop);
1156}
1157
1158////////////////////////////////////////////////////////////////////////////////
1159/// Exception safe assignement operator
1160
1162{
1163 TFormLeafInfoClones tmp(orig);
1164 Swap(tmp);
1165 return *this;
1166}
1167
1168////////////////////////////////////////////////////////////////////////////////
1169/// Return the current size of the the TClonesArray
1170
1172{
1173 if (!fCounter) {
1174 TClass *clonesClass = TClonesArray::Class();
1175 Int_t c_offset = 0;
1176 TStreamerElement *counter = ((TStreamerInfo*)clonesClass->GetStreamerInfo())->GetStreamerElement("fLast",c_offset);
1177 fCounter = new TFormLeafInfo(clonesClass,c_offset,counter);
1178 }
1179 return (Int_t)fCounter->ReadValue((char*)GetLocalValuePointer(leaf)) + 1;
1180}
1181
1182////////////////////////////////////////////////////////////////////////////////
1183/// Return the current size of the the TClonesArray
1184
1186{
1187 if (!fCounter) {
1188 TClass *clonesClass = TClonesArray::Class();
1189 Int_t c_offset = 0;
1190 TStreamerElement *counter = ((TStreamerInfo*)clonesClass->GetStreamerInfo())->GetStreamerElement("fLast",c_offset);
1191 fCounter = new TFormLeafInfo(clonesClass,c_offset,counter);
1192 }
1193 return (Int_t)fCounter->ReadValue(where) + 1;
1194}
1195
1196////////////////////////////////////////////////////////////////////////////////
1197/// Return the value of the underlying data member inside the
1198/// clones array.
1199
1200template <typename T>
1201T TFormLeafInfoClones::ReadValueImpl(char *where, Int_t instance)
1202{
1203 if (fNext==0) return 0;
1204 Int_t len,index,sub_instance;
1205 len = fNext->GetArrayLength();
1206 if (len) {
1207 index = instance / len;
1208 sub_instance = instance % len;
1209 } else {
1210 index = instance;
1211 sub_instance = 0;
1212 }
1213 TClonesArray * clones = (TClonesArray*)where;
1214 if (!clones) return 0;
1215 // Note we take advantage of having only one physically variable
1216 // dimension:
1217 char * obj = (char*)clones->UncheckedAt(index);
1218 return fNext->ReadTypedValue<T>(obj,sub_instance);
1219}
1220
1223
1224////////////////////////////////////////////////////////////////////////////////
1225/// Return the pointer to the clonesArray
1226
1228{
1229 TClonesArray * clones;
1230 if (fTop) {
1231 if (leaf->InheritsFrom(TLeafObject::Class()) ) {
1232 clones = (TClonesArray*)((TLeafObject*)leaf)->GetObject();
1233 } else {
1234 clones = (TClonesArray*)((TBranchElement*)leaf->GetBranch())->GetObject();
1235 }
1236 } else {
1238 }
1239 return clones;
1240}
1241
1242////////////////////////////////////////////////////////////////////////////////
1243/// Return the address of the underlying current value
1244
1246{
1248}
1249
1250////////////////////////////////////////////////////////////////////////////////
1251/// Return the value of the underlying data member inside the
1252/// clones array.
1253
1254template <typename T>
1255T TFormLeafInfoClones::GetValueImpl(TLeaf *leaf, Int_t instance)
1256{
1257 if (fNext==0) return 0;
1258 Int_t len,index,sub_instance;
1259 len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1260 Int_t primary = fNext->GetPrimaryIndex();
1261 if (len) {
1262 index = instance / len;
1263 sub_instance = instance % len;
1264 } else if (primary>=0) {
1265 index = primary;
1266 sub_instance = instance;
1267 } else {
1268 index = instance;
1269 sub_instance = 0;
1270 }
1272 if (clones==0) return 0;
1273
1274 // Note we take advantage of having only one physically variable
1275 // dimension:
1276 char * obj = (char*)clones->UncheckedAt(index);
1277 return fNext->ReadTypedValue<T>(obj,sub_instance);
1278}
1279
1280////////////////////////////////////////////////////////////////////////////////
1281/// Return the pointer to the clonesArray
1282
1284{
1286 if (fNext && clones) {
1287 // Same as in TFormLeafInfoClones::GetValue
1288 Int_t len,index,sub_instance;
1289 len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1290 if (len) {
1291 index = instance / len;
1292 sub_instance = instance % len;
1293 } else {
1294 index = instance;
1295 sub_instance = 0;
1296 }
1297 return fNext->GetValuePointer((char*)clones->UncheckedAt(index),
1298 sub_instance);
1299 }
1300 return clones;
1301}
1302
1303////////////////////////////////////////////////////////////////////////////////
1304/// Return the pointer to the clonesArray
1305
1307{
1308 TClonesArray * clones = (TClonesArray*) where;
1309 if (fNext) {
1310 // Same as in TFormLeafInfoClones::GetValue
1311 Int_t len,index,sub_instance;
1312 len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1313 if (len) {
1314 index = instance / len;
1315 sub_instance = instance % len;
1316 } else {
1317 index = instance;
1318 sub_instance = 0;
1319 }
1320 return fNext->GetValuePointer((char*)clones->UncheckedAt(index),
1321 sub_instance);
1322 }
1323 return clones;
1324}
1325
1326////////////////////////////////////////////////////////////////////////////////
1327/// \class TFormLeafInfoCollectionObject
1328/// A small helper class to implement reading a data member
1329/// on a TClonesArray object stored in a TTree.
1330
1331////////////////////////////////////////////////////////////////////////////////
1332/// Constructor.
1333
1335 TFormLeafInfo(classptr,0,R__GetFakeClonesElem()),fTop(top)
1336{
1337}
1338
1339////////////////////////////////////////////////////////////////////////////////
1340/// Constructor.
1341
1343 TFormLeafInfo(orig),fTop(orig.fTop)
1344{
1345}
1346
1347////////////////////////////////////////////////////////////////////////////////
1348/// Exception safe swap.
1349
1351{
1352 TFormLeafInfo::Swap(other);
1353 std::swap(fTop,other.fTop);
1354}
1355
1356////////////////////////////////////////////////////////////////////////////////
1357/// Exception safe assignement operator
1358
1360{
1362 Swap(tmp);
1363 return *this;
1364}
1365
1366////////////////////////////////////////////////////////////////////////////////
1367/// Return the current size of the the TClonesArray
1368
1370{
1371 return 1;
1372}
1373
1374////////////////////////////////////////////////////////////////////////////////
1375/// Return the value of the underlying data member inside the
1376/// clones array.
1377
1379{
1380 R__ASSERT(0);
1381 return 0;
1382}
1383
1384////////////////////////////////////////////////////////////////////////////////
1385/// Return the pointer to the clonesArray
1386
1388{
1389 void* collection;
1390 if (fTop) {
1391 if (leaf->InheritsFrom(TLeafObject::Class()) ) {
1392 collection = ((TLeafObject*)leaf)->GetObject();
1393 } else {
1394 collection = ((TBranchElement*)leaf->GetBranch())->GetObject();
1395 }
1396 } else {
1397 collection = TFormLeafInfo::GetLocalValuePointer(leaf);
1398 }
1399 return collection;
1400}
1401
1402////////////////////////////////////////////////////////////////////////////////
1403/// Return the address of the underlying current value
1404
1406{
1408}
1409
1410////////////////////////////////////////////////////////////////////////////////
1411/// Return the value of the underlying data member inside the
1412/// clones array.
1413
1414template <typename T>
1415T TFormLeafInfoCollectionObject::GetValueImpl(TLeaf *leaf, Int_t instance)
1416{
1417 char * obj = (char*)GetLocalValuePointer(leaf);
1418
1419 if (fNext==0) return 0;
1420 return fNext->ReadTypedValue<T>(obj,instance);
1421}
1422
1424
1425////////////////////////////////////////////////////////////////////////////////
1426/// Return the pointer to the clonesArray
1427
1429{
1430 void *collection = GetLocalValuePointer(leaf);
1431 if (fNext) {
1432 return fNext->GetValuePointer((char*)collection,instance);
1433 }
1434 return collection;
1435}
1436
1437////////////////////////////////////////////////////////////////////////////////
1438/// Return the pointer to the clonesArray
1439
1441{
1442 if (fNext) {
1443 return fNext->GetValuePointer(where,instance);
1444 }
1445 return where;
1446}
1447
1448////////////////////////////////////////////////////////////////////////////////
1449/// \class TFormLeafInfoCollection
1450/// A small helper class to implement reading a data
1451/// member on a generic collection object stored in a TTree.
1452
1453////////////////////////////////////////////////////////////////////////////////
1454/// Cosntructor.
1455
1457 Long_t offset,
1458 TStreamerElement* element,
1459 Bool_t top) :
1460 TFormLeafInfo(classptr,offset,element),
1461 fTop(top),
1462 fCollClass( 0),
1463 fCollProxy( 0),
1464 fLocalElement( 0)
1465{
1466 if (element) {
1467 fCollClass = element->GetClass();
1468 } else if (classptr) {
1469 fCollClass = classptr;
1470 }
1471 if (fCollClass
1474
1477 }
1478}
1479
1480////////////////////////////////////////////////////////////////////////////////
1481/// Constructor.
1482
1484 Long_t offset,
1485 TClass* elementclassptr,
1486 Bool_t top) :
1487 TFormLeafInfo(motherclassptr,offset,
1488 new TStreamerElement("collection","in class",
1489 0,
1490 TStreamerInfo::kAny,
1491 elementclassptr
1492 ? elementclassptr->GetName()
1493 : ( motherclassptr
1494 ? motherclassptr->GetName()
1495 : "Unknwon")
1496 ) ),
1497 fTop(top),
1498 fCollClass( 0),
1499 fCollProxy( 0) ,
1500 fLocalElement( fElement )
1501{
1502 if (elementclassptr) {
1503 fCollClass = elementclassptr;
1504 } else if (motherclassptr) {
1505 fCollClass = motherclassptr;
1506 }
1507 if (fCollClass
1510 {
1513 }
1514}
1515
1516////////////////////////////////////////////////////////////////////////////////
1517/// Constructor.
1518
1520 TFormLeafInfo(),
1521 fTop(kFALSE),
1522 fCollClass( 0),
1523 fCollProxy( 0),
1524 fLocalElement( 0)
1525{
1526}
1527
1528////////////////////////////////////////////////////////////////////////////////
1529/// Constructor.
1530
1532 TFormLeafInfo(orig),
1533 fTop( orig.fTop),
1534 fCollClass( orig.fCollClass ),
1535 fCollClassName( orig.fCollClassName ),
1536 fCollProxy( orig.fCollProxy ? orig.fCollProxy->Generate() : 0 ),
1537 fLocalElement( 0 ) // humm why not initialize it?
1538{
1539}
1540
1541////////////////////////////////////////////////////////////////////////////////
1542/// Exception safe swap.
1543
1545{
1546 TFormLeafInfo::Swap(other);
1547 std::swap(fTop,other.fTop);
1548 std::swap(fClass,other.fClass);
1552}
1553
1554////////////////////////////////////////////////////////////////////////////////
1555/// Exception safe assignment operator.
1556
1558{
1559 TFormLeafInfoCollection tmp(other);
1560 Swap(tmp);
1561 return *this;
1562}
1563
1564////////////////////////////////////////////////////////////////////////////////
1565/// Destructor.
1566
1568{
1569 delete fCollProxy;
1570 delete fLocalElement;
1571}
1572
1573////////////////////////////////////////////////////////////////////////////////
1574/// Copy of the object and its content.
1575
1577{
1578 return new TFormLeafInfoCollection(*this);
1579}
1580
1581////////////////////////////////////////////////////////////////////////////////
1582/// We reloading all cached information in case the underlying class
1583/// information has changed (for example when changing from the 'emulated'
1584/// class to the real class.
1585
1587{
1588 Bool_t changed = kFALSE;
1589 TClass * new_class = TClass::GetClass(fCollClassName);
1590 if (new_class!=fCollClass) {
1591 delete fCollProxy; fCollProxy = 0;
1592 fCollClass = new_class;
1595 }
1596 changed = kTRUE;
1597 }
1598 return changed || TFormLeafInfo::Update();
1599}
1600
1601////////////////////////////////////////////////////////////////////////////////
1602/// Return true if the underlying data has a array size counter
1603
1605{
1606 return fCounter!=0 || fCollProxy!=0;
1607}
1608
1609////////////////////////////////////////////////////////////////////////////////
1610/// Return the current size of the the TClonesArray
1611
1613{
1614 void *ptr = GetLocalValuePointer(leaf);
1615
1616 if (fCounter) { return (Int_t)fCounter->ReadValue((char*)ptr); }
1617
1619 if (ptr==0) return 0;
1621 return (Int_t)fCollProxy->Size();
1622}
1623
1624////////////////////////////////////////////////////////////////////////////////
1625/// Return the size of the underlying array for the current entry in the TTree.
1626
1628{
1629 if (fCounter) { return (Int_t)fCounter->ReadValue(where); }
1631 if (where==0) return 0;
1632 void *ptr = GetLocalValuePointer(where,0);
1634 return (Int_t)fCollProxy->Size();
1635}
1636
1637////////////////////////////////////////////////////////////////////////////////
1638/// Return the current size of the the TClonesArray
1639
1641{
1642 void *ptr = GetLocalValuePointer(leaf,instance);
1643 if (fCounter) {
1644 return (Int_t)fCounter->ReadValue((char*)ptr);
1645 }
1647 if (ptr==0) return 0;
1649 return (Int_t)fCollProxy->Size();
1650}
1651
1652////////////////////////////////////////////////////////////////////////////////
1653/// Return the value of the underlying data member inside the
1654/// clones array.
1655
1656template <typename T>
1657T TFormLeafInfoCollection::ReadValueImpl(char *where, Int_t instance)
1658{
1659 if (fNext==0) return 0;
1660 UInt_t len,index,sub_instance;
1661 len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1662 Int_t primary = fNext->GetPrimaryIndex();
1663 if (len) {
1664 index = instance / len;
1665 sub_instance = instance % len;
1666 } else if (primary>=0) {
1667 index = primary;
1668 sub_instance = instance;
1669 } else {
1670 index = instance;
1671 sub_instance = 0;
1672 }
1673
1675 void *ptr = GetLocalValuePointer(where,instance);
1677
1678 // Note we take advantage of having only one physically variable
1679 // dimension:
1680
1681 char * obj = (char*)fCollProxy->At(index);
1682 if (fCollProxy->HasPointers()) obj = *(char**)obj;
1683 return fNext->ReadTypedValue<T>(obj,sub_instance);
1684}
1685
1688
1689////////////////////////////////////////////////////////////////////////////////
1690/// Return the pointer to the clonesArray
1691
1693{
1694 void *collection;
1695 if (fTop) {
1696 if (leaf->InheritsFrom(TLeafObject::Class()) ) {
1697 collection = ((TLeafObject*)leaf)->GetObject();
1698 } else {
1699 collection = ((TBranchElement*)leaf->GetBranch())->GetObject();
1700 }
1701 } else {
1702 collection = TFormLeafInfo::GetLocalValuePointer(leaf);
1703 }
1704 return collection;
1705}
1706
1707////////////////////////////////////////////////////////////////////////////////
1708/// Return the address of the local value
1709
1711{
1713}
1714
1715////////////////////////////////////////////////////////////////////////////////
1716/// Return the value of the underlying data member inside the
1717/// clones array.
1718
1719template <typename T>
1720T TFormLeafInfoCollection::GetValueImpl(TLeaf *leaf, Int_t instance)
1721{
1722 if (fNext==0) return 0;
1723 Int_t len,index,sub_instance;
1724 len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1725 Int_t primary = fNext->GetPrimaryIndex();
1726 if (len) {
1727 index = instance / len;
1728 sub_instance = instance % len;
1729 } else if (primary>=0) {
1730 index = primary;
1731 sub_instance = instance;
1732 } else {
1733 index = instance;
1734 sub_instance = 0;
1735 }
1736
1738 void *coll = GetLocalValuePointer(leaf);
1740
1741 // Note we take advantage of having only one physically variable
1742 // dimension:
1743 char * obj = (char*)fCollProxy->At(index);
1744 if (obj==0) return 0;
1745 if (fCollProxy->HasPointers()) obj = *(char**)obj;
1746 if (obj==0) return 0;
1747 return fNext->ReadTypedValue<T>(obj,sub_instance);
1748}
1749
1750////////////////////////////////////////////////////////////////////////////////
1751/// Return the pointer to the clonesArray
1752
1754{
1756
1757 void *collection = GetLocalValuePointer(leaf);
1758
1759 if (fNext) {
1760 // Same as in TFormLeafInfoClones::GetValue
1761 Int_t len,index,sub_instance;
1762 if (fNext->fElement &&
1763 (fNext->fNext || !fNext->IsString()) ) {
1764 len = fNext->GetArrayLength();
1765 } else {
1766 len = 0;
1767 }
1768 if (len) {
1769 index = instance / len;
1770 sub_instance = instance % len;
1771 } else {
1772 index = instance;
1773 sub_instance = 0;
1774 }
1776 char * obj = (char*)fCollProxy->At(index);
1777 if (fCollProxy->HasPointers()) obj = *(char**)obj;
1778 return fNext->GetValuePointer(obj,sub_instance);
1779 }
1780 return collection;
1781}
1782
1783////////////////////////////////////////////////////////////////////////////////
1784/// Return the pointer to the clonesArray
1785
1787{
1789
1790 void *collection = where;
1791
1792 if (fNext) {
1793 // Same as in TFormLeafInfoClones::GetValue
1794 Int_t len,index,sub_instance;
1795 len = (fNext->fElement==0)? 0 : fNext->GetArrayLength();
1796 if (len) {
1797 index = instance / len;
1798 sub_instance = instance % len;
1799 } else {
1800 index = instance;
1801 sub_instance = 0;
1802 }
1804 char * obj = (char*)fCollProxy->At(index);
1805 if (fCollProxy->HasPointers()) obj = *(char**)obj;
1806 return fNext->GetValuePointer(obj,sub_instance);
1807 }
1808 return collection;
1809}
1810
1811////////////////////////////////////////////////////////////////////////////////
1812/// \class TFormLeafInfoCollectionSize
1813/// Used to return the size of a collection
1814
1815////////////////////////////////////////////////////////////////////////////////
1816/// Constructor.
1817
1819 TFormLeafInfo(), fCollClass(classptr), fCollProxy(0)
1820{
1821 if (fCollClass
1824
1827 }
1828}
1829
1830////////////////////////////////////////////////////////////////////////////////
1831/// Constructor.
1832
1834 TClass* classptr,Long_t offset,TStreamerElement* element) :
1835 TFormLeafInfo(classptr,offset,element), fCollClass(element->GetClassPointer()), fCollProxy(0)
1836{
1837 if (fCollClass
1840
1843 }
1844}
1845
1846////////////////////////////////////////////////////////////////////////////////
1847/// Constructor.
1848
1850 TFormLeafInfo(), fCollClass(0), fCollProxy(0)
1851{
1852}
1853
1854////////////////////////////////////////////////////////////////////////////////
1855/// Constructor.
1856
1859 fCollClass(orig.fCollClass),
1860 fCollClassName(orig.fCollClassName),
1861 fCollProxy(orig.fCollProxy?orig.fCollProxy->Generate():0)
1862{
1863}
1864
1865////////////////////////////////////////////////////////////////////////////////
1866/// Exception safe assignment operator.
1867
1869{
1870 TFormLeafInfoCollectionSize tmp(other);
1871 Swap(tmp);
1872 return *this;
1873}
1874
1875////////////////////////////////////////////////////////////////////////////////
1876/// Exception safe swap.
1877
1879{
1880 TFormLeafInfo::Swap(other);
1884}
1885
1886////////////////////////////////////////////////////////////////////////////////
1887/// Destructor.
1888
1890{
1891 delete fCollProxy;
1892}
1893
1894////////////////////////////////////////////////////////////////////////////////
1895/// Copy the object and all of its content.
1896
1898{
1899 return new TFormLeafInfoCollectionSize(*this);
1900}
1901
1902////////////////////////////////////////////////////////////////////////////////
1903/// We reloading all cached information in case the underlying class
1904/// information has changed (for example when changing from the 'emulated'
1905/// class to the real class.
1906
1908{
1909 Bool_t changed = kFALSE;
1911 if (new_class!=fCollClass) {
1912 delete fCollProxy; fCollProxy = 0;
1913 fCollClass = new_class;
1916 }
1917 changed = kTRUE;
1918 }
1919 return changed;
1920}
1921
1922////////////////////////////////////////////////////////////////////////////////
1923/// Not implemented.
1924
1926{
1927 Error("GetValuePointer","This should never be called");
1928 return 0;
1929}
1930
1931////////////////////////////////////////////////////////////////////////////////
1932/// Not implemented.
1933
1934void *TFormLeafInfoCollectionSize::GetValuePointer(char * /* from */, Int_t /* instance */)
1935{
1936 Error("GetValuePointer","This should never be called");
1937 return 0;
1938}
1939
1940////////////////////////////////////////////////////////////////////////////////
1941/// Not implemented.
1942
1944{
1945 Error("GetLocalValuePointer","This should never be called");
1946 return 0;
1947}
1948
1949////////////////////////////////////////////////////////////////////////////////
1950/// Not implemented.
1951
1952void *TFormLeafInfoCollectionSize::GetLocalValuePointer( char * /* from */, Int_t /* instance */)
1953{
1954 Error("GetLocalValuePointer","This should never be called");
1955 return 0;
1956}
1957
1958////////////////////////////////////////////////////////////////////////////////
1959/// Return the value of the underlying pointer data member
1960
1962{
1964 if (where==0) return 0;
1965 void *ptr = fElement ? TFormLeafInfo::GetLocalValuePointer(where) : where;
1967 return (Int_t)fCollProxy->Size();
1968}
1969
1970////////////////////////////////////////////////////////////////////////////////
1971/// \class TFormLeafInfoPointer
1972/// A small helper class to implement reading a data
1973/// member by following a pointer inside a branch of TTree.
1974
1975////////////////////////////////////////////////////////////////////////////////
1976/// Constructor.
1977
1979 Long_t offset,
1980 TStreamerElement* element) :
1981 TFormLeafInfo(classptr,offset,element)
1982{
1983}
1984
1985////////////////////////////////////////////////////////////////////////////////
1986/// Copy the object and all of its contnet.
1987
1989{
1990 return new TFormLeafInfoPointer(*this);
1991}
1992
1993////////////////////////////////////////////////////////////////////////////////
1994/// Return the value of the underlying pointer data member
1995
1996template <typename T>
1997T TFormLeafInfoPointer::ReadValueImpl(char *where, Int_t instance)
1998{
1999 if (!fNext) return 0;
2000 char * whereoffset = where+fOffset;
2001 switch (fElement->GetNewType()) {
2002 // basic types
2008 {TObject **obj = (TObject**)(whereoffset);
2009 return obj && *obj ? fNext->ReadTypedValue<T>((char*)*obj,instance) : 0; }
2010
2018 {
2019 TObject *obj = (TObject*)(whereoffset);
2020 return fNext->ReadTypedValue<T>((char*)obj,instance);
2021 }
2022
2026 {
2027 Int_t len, index, sub_instance;
2028
2029 if (fNext) len = fNext->GetArrayLength();
2030 else len = 1;
2031 if (len) {
2032 index = instance / len;
2033 sub_instance = instance % len;
2034 } else {
2035 index = instance;
2036 sub_instance = 0;
2037 }
2038
2039 whereoffset += index*fElement->GetClassPointer()->Size();
2040
2041 TObject *obj = (TObject*)(whereoffset);
2042 return fNext->ReadTypedValue<T>((char*)obj,sub_instance);
2043 }
2044
2050 {
2051 TObject *obj = (TObject*)(whereoffset);
2052 return fNext->ReadTypedValue<T>((char*)obj,instance);
2053 }
2054
2055 case kOther_t:
2056 default: return 0;
2057 }
2058
2059}
2060
2061////////////////////////////////////////////////////////////////////////////////
2062/// Return the value of the underlying pointer data member
2063
2064template <typename T>
2065T TFormLeafInfoPointer::GetValueImpl(TLeaf *leaf, Int_t instance)
2066{
2067 if (!fNext) return 0;
2068 char * where = (char*)GetLocalValuePointer(leaf,instance);
2069 if (where==0) return 0;
2070 return fNext->ReadTypedValue<T>(where,instance);
2071}
2072
2075
2076////////////////////////////////////////////////////////////////////////////////
2077/// \class TFormLeafInfoMethod
2078/// Asmall helper class to implement executing a method
2079/// of an object stored in a TTree
2080
2081////////////////////////////////////////////////////////////////////////////////
2082/// Constructor.
2083
2085 TMethodCall *method) :
2086 TFormLeafInfo(classptr,0,0),fMethod(method),
2087 fResult(0), fCopyFormat(),fDeleteFormat(),fValuePointer(0),fIsByValue(kFALSE)
2088{
2089 if (method) {
2090 fMethodName = method->GetMethodName();
2091 fParams = method->GetParams();
2093 if (r == TMethodCall::kOther) {
2094 const char* rtype = fMethod->GetMethod()->GetReturnTypeName();
2095 Long_t rprop = fMethod->GetMethod()->Property();
2096 if (rtype[strlen(rtype)-1]!='*' &&
2097 rtype[strlen(rtype)-1]!='&' &&
2098 !(rprop & (kIsPointer|kIsReference)) ) {
2099 fCopyFormat = "new ";
2100 fCopyFormat += rtype;
2101 fCopyFormat += "(*(";
2102 fCopyFormat += rtype;
2103 fCopyFormat += "*)0x%lx)";
2104
2105 fDeleteFormat = "delete (";
2106 fDeleteFormat += rtype;
2107 fDeleteFormat += "*)0x%lx";
2108
2109 fIsByValue = kTRUE;
2110 }
2111 }
2112 }
2113}
2114
2115////////////////////////////////////////////////////////////////////////////////
2116/// Constructor.
2117
2119 : TFormLeafInfo(orig)
2120{
2121 fMethodName = orig.fMethodName;
2122 fParams = orig.fParams ;
2123 fResult = orig.fResult;
2124 if (orig.fMethod) {
2125 fMethod = new TMethodCall();
2126 fMethod->Init(orig.fMethod->GetMethod());
2127 } else {
2128 fMethod = 0;
2129 }
2130 fCopyFormat = orig.fCopyFormat;
2132 fValuePointer = 0;
2133 fIsByValue = orig.fIsByValue;
2134}
2135
2136
2137////////////////////////////////////////////////////////////////////////////////
2138/// Exception safe swap.
2139
2141{
2142 TFormLeafInfo::Swap(other);
2143 std::swap(fMethod,other.fMethod);
2145 std::swap(fParams,other.fParams);
2146 std::swap(fResult,other.fResult);
2151}
2152
2153////////////////////////////////////////////////////////////////////////////////
2154/// Exception safe assignment operator.
2155
2157{
2158 TFormLeafInfoMethod tmp(other);
2159 Swap(tmp);
2160 return *this;
2161}
2162
2163////////////////////////////////////////////////////////////////////////////////
2164/// Destructor.
2165
2167{
2168 if (fValuePointer) {
2170 }
2171 delete fMethod;
2172}
2173
2174////////////////////////////////////////////////////////////////////////////////
2175/// Copy the object and all its content.
2176
2178{
2179 return new TFormLeafInfoMethod(*this);
2180}
2181
2182////////////////////////////////////////////////////////////////////////////////
2183/// Return the TClass corresponding to the return type of the function
2184/// if it is an object type or if the return type is a reference (&) then
2185/// return the TClass corresponding to the referencee.
2186
2188{
2189 if (!mc || !mc->GetMethod())
2190 return nullptr;
2191
2192 std::string return_type;
2193
2194 if (0 == strcmp(mc->GetMethod()->GetReturnTypeName(), "void"))
2195 return nullptr;
2196
2198
2199 {
2202 }
2203 // Beyhond this point we no longer 'need' the lock.
2204 // How TClass::GetClass will take at least the read lock to search
2205 // So keeping it just a little longer is likely to be faster
2206 // than releasing and retaking it.
2207
2208 return_type = gInterpreter->TypeName(return_type.c_str());
2209
2210 if (return_type == "void")
2211 return nullptr;
2212
2213 return TClass::GetClass(return_type.c_str());
2214}
2215
2216////////////////////////////////////////////////////////////////////////////////
2217/// Return the type of the underlying return value
2218
2220{
2221 if (fNext) return fNext->GetClass();
2223 if (r!=TMethodCall::kOther) return 0;
2224
2225 return ReturnTClass(fMethod);
2226}
2227
2228////////////////////////////////////////////////////////////////////////////////
2229/// Return true if the return value is integral.
2230
2232{
2234 if (r == TMethodCall::kLong) {
2235 return kTRUE;
2236 } else return kFALSE;
2237}
2238
2239////////////////////////////////////////////////////////////////////////////////
2240/// Return true if the return value is a string.
2241
2243{
2244 if (fNext) return fNext->IsString();
2245
2247 return (r==TMethodCall::kString);
2248}
2249
2250////////////////////////////////////////////////////////////////////////////////
2251/// We reloading all cached information in case the underlying class
2252/// information has changed (for example when changing from the 'emulated'
2253/// class to the real class.
2254
2256{
2257 if (!TFormLeafInfo::Update()) return kFALSE;
2258 delete fMethod;
2260 return kTRUE;
2261}
2262
2263////////////////////////////////////////////////////////////////////////////////
2264/// This is implemented here because some compiler want ALL the
2265/// signature of an overloaded function to be re-implemented.
2266
2269{
2271}
2272
2273////////////////////////////////////////////////////////////////////////////////
2274/// Return the address of the lcoal underlying value.
2275
2277 Int_t /*instance*/)
2278{
2279 void *thisobj = from;
2280 if (!thisobj) return 0;
2281
2283 fResult = 0;
2284
2285 if (r == TMethodCall::kLong) {
2286 Long_t l = 0;
2287 fMethod->Execute(thisobj, l);
2288 fResult = (Double_t) l;
2289 // Get rid of temporary return object.
2290 gInterpreter->ClearStack();
2291 return &fResult;
2292
2293 } else if (r == TMethodCall::kDouble) {
2294 Double_t d = 0;
2295 fMethod->Execute(thisobj, d);
2296 fResult = (Double_t) d;
2297 // Get rid of temporary return object.
2298 gInterpreter->ClearStack();
2299 return &fResult;
2300
2301 } else if (r == TMethodCall::kString) {
2302 char *returntext = 0;
2303 fMethod->Execute(thisobj,&returntext);
2304 gInterpreter->ClearStack();
2305 return returntext;
2306
2307 } else if (r == TMethodCall::kOther) {
2308 char * char_result = 0;
2309 if (fIsByValue) {
2310 if (fValuePointer) {
2311 gROOT->ProcessLine(Form(fDeleteFormat.Data(),fValuePointer));
2312 fValuePointer = 0;
2313 }
2314 }
2315 fMethod->Execute(thisobj, &char_result);
2316 if (fIsByValue) {
2317 fValuePointer = (char*)gInterpreter->Calc(Form(fCopyFormat.Data(),char_result));
2318 char_result = (char*)fValuePointer;
2319 }
2320 gInterpreter->ClearStack();
2321 return char_result;
2322
2323 }
2324 return 0;
2325}
2326
2327////////////////////////////////////////////////////////////////////////////////
2328/// Execute the method on the given address
2329
2330template <typename T>
2331T TFormLeafInfoMethod::ReadValueImpl(char *where, Int_t instance)
2332{
2333 void *thisobj = where;
2334 if (!thisobj) return 0;
2335
2337 T result = 0;
2338
2339 if (r == TMethodCall::kLong) {
2340 Long_t l = 0;
2341 fMethod->Execute(thisobj, l);
2342 result = (T) l;
2343
2344 } else if (r == TMethodCall::kDouble) {
2345 Double_t d = 0;
2346 fMethod->Execute(thisobj, d);
2347 result = (T) d;
2348
2349 } else if (r == TMethodCall::kString) {
2350 char *returntext = 0;
2351 fMethod->Execute(thisobj,&returntext);
2352 result = T((Long_t) returntext);
2353
2354 } else if (fNext) {
2355 char * char_result = 0;
2356 fMethod->Execute(thisobj, &char_result);
2357 result = fNext->ReadTypedValue<T>(char_result,instance);
2358
2359 } else fMethod->Execute(thisobj);
2360
2361 // Get rid of temporary return object.
2362 gInterpreter->ClearStack();
2363 return result;
2364}
2365
2367
2368////////////////////////////////////////////////////////////////////////////////
2369/// \class TFormLeafInfoMultiVarDim
2370/// A helper class to implement reading a
2371/// data member on a variable size array inside a TClonesArray object stored in
2372/// a TTree. This is the version used when the data member is inside a
2373/// non-split object.
2374
2375////////////////////////////////////////////////////////////////////////////////
2376/// Constructor.
2377
2379 Long_t offset,
2380 TStreamerElement* element,
2381 TFormLeafInfo* parent) :
2382 TFormLeafInfo(classptr,offset,element),fNsize(0),fCounter2(0),fSumOfSizes(0),
2383 fDim(0),fVirtDim(-1),fPrimaryIndex(-1),fSecondaryIndex(-1)
2384{
2385 if (element && element->InheritsFrom(TStreamerBasicPointer::Class())) {
2387
2388 Int_t counterOffset = 0;
2389 TStreamerElement* counter = ((TStreamerInfo*)classptr->GetStreamerInfo())->GetStreamerElement(elem->GetCountName(),counterOffset);
2390 if (!parent) return;
2391 fCounter2 = parent->DeepCopy();
2392 TFormLeafInfo ** next = &(fCounter2->fNext);
2393 while(*next != 0) next = &( (*next)->fNext);
2394 *next = new TFormLeafInfo(classptr,counterOffset,counter);
2395
2396 } else Error("Constructor","Called without a proper TStreamerElement");
2397}
2398
2399////////////////////////////////////////////////////////////////////////////////
2400/// Constructor.
2401
2403 TFormLeafInfo(0,0,0),fNsize(0),fCounter2(0),fSumOfSizes(0),
2404 fDim(0),fVirtDim(-1),fPrimaryIndex(-1),fSecondaryIndex(-1)
2405{
2406}
2407
2408////////////////////////////////////////////////////////////////////////////////
2409/// Constructor.
2410
2412{
2413 fNsize = orig.fNsize;
2414 orig.fSizes.Copy(fSizes);
2415 fCounter2 = orig.fCounter2?orig.fCounter2->DeepCopy():0;
2416 fSumOfSizes = orig.fSumOfSizes;
2417 fDim = orig.fDim;
2418 fVirtDim = orig.fVirtDim;
2421}
2422
2423////////////////////////////////////////////////////////////////////////////////
2424/// Exception safe swap.
2425
2427{
2428 TFormLeafInfo::Swap(other);
2429 std::swap(fNsize,other.fNsize);
2430 std::swap(fSizes,other.fSizes);
2432 std::swap(fDim,other.fDim);
2436}
2437
2438////////////////////////////////////////////////////////////////////////////////
2439/// Exception safe assignment operator.
2440
2442{
2443 TFormLeafInfoMultiVarDim tmp(other);
2444 Swap(tmp);
2445 return *this;
2446}
2447
2448////////////////////////////////////////////////////////////////////////////////
2449/// Copy the object and all its content.
2450
2452{
2453 return new TFormLeafInfoMultiVarDim(*this);
2454}
2455
2456////////////////////////////////////////////////////////////////////////////////
2457/// Destructor.
2458
2460{
2461 delete fCounter2;
2462}
2463
2464/* The proper indexing and unwinding of index is done by prior leafinfo in the chain. */
2465//virtual Double_t TFormLeafInfoMultiVarDim::ReadValue(char *where, Int_t instance = 0) {
2466// return TFormLeafInfo::ReadValue(where,instance);
2467//}
2468
2469////////////////////////////////////////////////////////////////////////////////
2470/// Load the current array sizes.
2471
2473{
2474 if (fElement) {
2475 TLeaf *leaf = (TLeaf*)branch->GetListOfLeaves()->At(0);
2476 if (fCounter) fNsize = (Int_t)fCounter->GetValue(leaf);
2477 else fNsize = fCounter2->GetCounterValue(leaf);
2479 fSumOfSizes = 0;
2480 for (Int_t i=0; i<fNsize; i++) {
2481 Int_t size = (Int_t)fCounter2->GetValue(leaf,i);
2482 fSumOfSizes += size;
2483 fSizes.AddAt( size, i );
2484 }
2485 return;
2486 }
2487 if (!fCounter2 || !fCounter) return;
2488 TBranchElement *br = dynamic_cast<TBranchElement*>(branch);
2489 R__ASSERT(br);
2490 fNsize = br->GetBranchCount()->GetNdata();
2492 fSumOfSizes = 0;
2493 for (Int_t i=0; i<fNsize; i++) {
2494 Int_t size = (Int_t)fCounter2->GetValue((TLeaf*)br->GetBranchCount2()->GetListOfLeaves()->At(0),i);
2495 fSumOfSizes += size;
2496 fSizes.AddAt( size, i );
2497 }
2498}
2499
2500////////////////////////////////////////////////////////////////////////////////
2501/// Return the index vlaue of the primary index.
2502
2504{
2505 return fPrimaryIndex;
2506}
2507
2508////////////////////////////////////////////////////////////////////////////////
2509/// Return the size of the requested sub-array.
2510
2512{
2513 if (index >= fSizes.GetSize()) {
2514 return -1;
2515 } else {
2516 return fSizes.At(index);
2517 }
2518}
2519
2520////////////////////////////////////////////////////////////////////////////////
2521/// Set the current value of the primary index.
2522
2524{
2525 fPrimaryIndex = index;
2526}
2527
2528////////////////////////////////////////////////////////////////////////////////
2529/// Set the current value of the primary index.
2530
2532{
2533 fSecondaryIndex = index;
2534}
2535
2536////////////////////////////////////////////////////////////////////////////////
2537/// Set the sizes of the sub-array.
2538
2540{
2541 fSumOfSizes += (val - fSizes.At(index));
2542 fSizes.AddAt(val,index);
2543}
2544
2545////////////////////////////////////////////////////////////////////////////////
2546/// Get the total size.
2547
2549{
2550 return fSumOfSizes;
2551}
2552
2553////////////////////////////////////////////////////////////////////////////////
2554
2556 Int_t /*instance*/)
2557{
2558 /* The proper indexing and unwinding of index need to be done by prior leafinfo in the chain. */
2559 Error("GetValue","This should never be called");
2560 return 0;
2561}
2562
2563
2564////////////////////////////////////////////////////////////////////////////////
2565/// Return the index of the dimension which varies
2566/// for each elements of an enclosing array (typically a TClonesArray)
2567
2569{
2570 return fDim;
2571}
2572
2573////////////////////////////////////////////////////////////////////////////////
2574/// Return the virtual index (for this expression) of the dimension which varies
2575/// for each elements of an enclosing array (typically a TClonesArray)
2576
2578{
2579 return fVirtDim;
2580}
2581
2582////////////////////////////////////////////////////////////////////////////////
2583/// We reloading all cached information in case the underlying class
2584/// information has changed (for example when changing from the 'emulated'
2585/// class to the real class.
2586
2588{
2590 if (fCounter2) fCounter2->Update();
2591 return res;
2592}
2593
2594////////////////////////////////////////////////////////////////////////////////
2595/// Update the sizes of the arrays.
2596
2598{
2599 if (!garr) return;
2600 if (garr->GetSize()<fNsize) garr->Set(fNsize);
2601 for (Int_t i=0; i<fNsize; i++) {
2602 Int_t local = fSizes.At(i);
2603 Int_t global = garr->At(i);
2604 if (global==0 || local<global) global = local;
2605 garr->AddAt(global,i);
2606 }
2607}
2608
2609////////////////////////////////////////////////////////////////////////////////
2610/// \class TFormLeafInfoMultiVarDimDirect
2611/// A small helper class to implement reading
2612/// a data member on a variable size array inside a TClonesArray object stored
2613/// in a TTree. This is the version used for split access
2614
2615////////////////////////////////////////////////////////////////////////////////
2616/// Copy the object and all its content.
2617
2619{
2620 return new TFormLeafInfoMultiVarDimDirect(*this);
2621}
2622
2623////////////////////////////////////////////////////////////////////////////////
2624/// Return the undersying value.
2625
2626template <typename T>
2627T TFormLeafInfoMultiVarDimDirect::GetValueImpl(TLeaf *leaf, Int_t instance)
2628{
2629 return ((TLeafElement*)leaf)->GetTypedValueSubArray<T>(fPrimaryIndex,instance);
2630}
2631
2633
2634////////////////////////////////////////////////////////////////////////////////
2635/// Not implemented.
2636
2638{
2639 Error("ReadValue","This should never be called");
2640 return 0;
2641}
2642
2643////////////////////////////////////////////////////////////////////////////////
2644/// \class TFormLeafInfoMultiVarDimCollection
2645/// A small helper class to implement reading
2646/// a data member on a variable size array inside a TClonesArray object stored
2647/// in a TTree. This is the version used for split access
2648
2649////////////////////////////////////////////////////////////////////////////////
2650/// Constructor.
2651
2653 TClass* motherclassptr,
2654 Long_t offset,
2655 TClass* elementclassptr,
2656 TFormLeafInfo *parent) :
2657 TFormLeafInfoMultiVarDim(motherclassptr,offset,
2658 new TStreamerElement("collection","in class",
2659 0,
2660 TStreamerInfo::kAny,
2661 elementclassptr
2662 ? elementclassptr->GetName()
2663 : ( motherclassptr
2664 ? motherclassptr->GetName()
2665 : "Unknwon")
2666 )
2667 )
2668{
2669 R__ASSERT(parent);
2670 fCounter = parent->DeepCopy();
2671 fCounter2 = parent->DeepCopy();
2672 TFormLeafInfo ** next = &(fCounter2->fNext);
2673 while(*next != 0) next = &( (*next)->fNext);
2674 *next = new TFormLeafInfoCollectionSize(elementclassptr);
2675}
2676
2677////////////////////////////////////////////////////////////////////////////////
2678/// Constructor.
2679
2681 TClass* motherclassptr,
2682 Long_t offset,
2683 TStreamerElement* element,
2684 TFormLeafInfo *parent) :
2685 TFormLeafInfoMultiVarDim(motherclassptr,offset,element)
2686{
2687 R__ASSERT(parent && element);
2688 fCounter = parent->DeepCopy();
2689 fCounter2 = parent->DeepCopy();
2690 TFormLeafInfo ** next = &(fCounter2->fNext);
2691 while(*next != 0) next = &( (*next)->fNext);
2692 *next = new TFormLeafInfoCollectionSize(motherclassptr,offset,element);
2693}
2694
2695////////////////////////////////////////////////////////////////////////////////
2696/// Copy the object and all its content.
2697
2699{
2700 return new TFormLeafInfoMultiVarDimCollection(*this);
2701}
2702
2703////////////////////////////////////////////////////////////////////////////////
2704
2706 Int_t /* instance */)
2707{
2708 /* The proper indexing and unwinding of index need to be done by prior leafinfo in the chain. */
2709 Error("GetValue","This should never be called");
2710 return 0;
2711}
2712
2713////////////////////////////////////////////////////////////////////////////////
2714/// Load the current array sizes.
2715
2717{
2719
2720 TLeaf *leaf = (TLeaf*)branch->GetListOfLeaves()->At(0);
2722
2724 fSumOfSizes = 0;
2725 for (Int_t i=0; i<fNsize; i++) {
2726 Int_t size = (Int_t)fCounter2->GetValue(leaf,i);
2727 fSumOfSizes += size;
2728 fSizes.AddAt( size, i );
2729 }
2730 return;
2731}
2732
2733////////////////////////////////////////////////////////////////////////////////
2734/// Return the value of the underlying data.
2735
2736template <typename T>
2737T TFormLeafInfoMultiVarDimCollection::ReadValueImpl(char *where, Int_t instance)
2738{
2739 if (fSecondaryIndex>=0) {
2740 UInt_t len = fNext->GetArrayLength();
2741 if (len) {
2743 } else {
2745 }
2746 }
2747 return fNext->ReadTypedValue<T>(where,instance);
2748}
2749
2751
2752
2753////////////////////////////////////////////////////////////////////////////////
2754/// \class TFormLeafInfoMultiVarDimClones
2755/// A small helper class to implement reading
2756/// a data member on a variable size array inside a TClonesArray object stored
2757/// in a TTree. This is the version used for split access
2758
2759////////////////////////////////////////////////////////////////////////////////
2760/// Constructor.
2761
2763 TClass* motherclassptr,
2764 Long_t offset,
2765 TClass* elementclassptr,
2766 TFormLeafInfo *parent) :
2767 TFormLeafInfoMultiVarDim(motherclassptr,offset,
2768 new TStreamerElement("clones","in class",
2769 0,
2770 TStreamerInfo::kAny,
2771 elementclassptr
2772 ? elementclassptr->GetName()
2773 : ( motherclassptr
2774 ? motherclassptr->GetName()
2775 : "Unknwon")
2776 )
2777 )
2778{
2779 R__ASSERT(parent);
2780 fCounter = parent->DeepCopy();
2781 fCounter2 = parent->DeepCopy();
2782 TFormLeafInfo ** next = &(fCounter2->fNext);
2783 while(*next != 0) next = &( (*next)->fNext);
2784 *next = new TFormLeafInfoClones(elementclassptr);
2785}
2786
2787////////////////////////////////////////////////////////////////////////////////
2788/// Constructor.
2789
2791 TClass* motherclassptr,
2792 Long_t offset,
2793 TStreamerElement* element,
2794 TFormLeafInfo *parent) :
2795 TFormLeafInfoMultiVarDim(motherclassptr,offset,element)
2796{
2797 R__ASSERT(parent && element);
2798 fCounter = parent->DeepCopy();
2799 fCounter2 = parent->DeepCopy();
2800 TFormLeafInfo ** next = &(fCounter2->fNext);
2801 while(*next != 0) next = &( (*next)->fNext);
2802 *next = new TFormLeafInfoClones(motherclassptr,offset,element);
2803}
2804
2805////////////////////////////////////////////////////////////////////////////////
2806/// Copy the object and all its data.
2807
2809{
2810 return new TFormLeafInfoMultiVarDimClones(*this);
2811}
2812
2813////////////////////////////////////////////////////////////////////////////////
2814
2816 Int_t /* instance */)
2817{
2818 /* The proper indexing and unwinding of index need to be done by prior leafinfo in the chain. */
2819 Error("GetValue","This should never be called");
2820 return 0;
2821}
2822
2823////////////////////////////////////////////////////////////////////////////////
2824/// Load the current array sizes.
2825
2827{
2829
2830 TLeaf *leaf = (TLeaf*)branch->GetListOfLeaves()->At(0);
2832
2834 fSumOfSizes = 0;
2835 for (Int_t i=0; i<fNsize; i++) {
2837 if (clones) {
2838 Int_t size = clones->GetEntries();
2839 fSumOfSizes += size;
2840 fSizes.AddAt( size, i );
2841 }
2842 }
2843 return;
2844}
2845
2846////////////////////////////////////////////////////////////////////////////////
2847/// Return the value of the underlying data.
2848
2849template <typename T>
2850T TFormLeafInfoMultiVarDimClones::ReadValueImpl(char *where, Int_t instance)
2851{
2852 if (fSecondaryIndex>=0) {
2853 UInt_t len = fNext->GetArrayLength();
2854 if (len) {
2856 } else {
2858 }
2859 }
2860 return fNext->ReadTypedValue<T>(where,instance);
2861}
2862
2864
2865////////////////////////////////////////////////////////////////////////////////
2866/// \class TFormLeafInfoCast
2867/// A small helper class to implement casting an object to
2868/// a different type (equivalent to dynamic_cast)
2869
2870////////////////////////////////////////////////////////////////////////////////
2871/// Constructor.
2872
2874 TFormLeafInfo(classptr),fCasted(casted),fGoodCast(kTRUE)
2875{
2876 if (casted) { fCastedName = casted->GetName(); }
2877 fMultiplicity = -1;
2879}
2880
2881////////////////////////////////////////////////////////////////////////////////
2882/// Constructor.
2883
2885 TFormLeafInfo(orig)
2886{
2887 fCasted = orig.fCasted;
2888 fCastedName = orig.fCastedName;
2889 fGoodCast = orig.fGoodCast;
2890 fIsTObject = orig.fIsTObject;
2891}
2892
2893////////////////////////////////////////////////////////////////////////////////
2894/// Exception safe swap.
2895
2897{
2898 TFormLeafInfo::Swap(other);
2899 std::swap(fCasted,other.fCasted);
2903}
2904
2905////////////////////////////////////////////////////////////////////////////////
2906/// Exception safe assignment operator.
2907
2909{
2910 TFormLeafInfoCast tmp(other);
2911 Swap(tmp);
2912 return *this;
2913}
2914
2915////////////////////////////////////////////////////////////////////////////////
2916/// Copy the object and all its content.
2917
2919{
2920 return new TFormLeafInfoCast(*this);
2921}
2922
2923////////////////////////////////////////////////////////////////////////////////
2924/// Destructor.
2925
2927{
2928}
2929
2930// Currently only implemented in TFormLeafInfoCast
2932{
2933 // Get the number of element in the entry.
2934
2935 if (!fGoodCast) return 0;
2936 if (fNext) return fNext->GetNdata();
2937 return 1;
2938}
2939
2940////////////////////////////////////////////////////////////////////////////////
2941/// Read the value at the given memory location
2942
2943template <typename T>
2944T TFormLeafInfoCast::ReadValueImpl(char *where, Int_t instance)
2945{
2946 if (!fNext) return 0;
2947
2948 // First check that the real class inherits from the
2949 // casted class
2950 // First assume TObject ...
2951 if ( fIsTObject && !((TObject*)where)->InheritsFrom(fCasted) ) {
2952 fGoodCast = kFALSE;
2953 return 0;
2954 } else {
2955 // We know we have a TBranchElement and we need to find out the
2956 // real class name.
2957 }
2958 fGoodCast = kTRUE;
2959 return fNext->ReadTypedValue<T>(where,instance);
2960}
2961
2963
2964
2965////////////////////////////////////////////////////////////////////////////////
2966/// We reloading all cached information in case the underlying class
2967/// information has changed (for example when changing from the 'emulated'
2968/// class to the real class.
2969
2971{
2972 if (fCasted) {
2973 TClass * new_class = TClass::GetClass(fCastedName);
2974 if (new_class!=fCasted) {
2975 fCasted = new_class;
2976 }
2977 }
2978 return TFormLeafInfo::Update();
2979}
2980
2981////////////////////////////////////////////////////////////////////////////////
2982/// \class TFormLeafInfoTTree
2983/// A small helper class to implement reading
2984/// from the containing TTree object itself.
2985
2987TFormLeafInfo( TTree::Class(), 0, 0 ), fTree(tree),fCurrent(current),fAlias(alias)
2988{
2989 // Constructor.
2990
2991 if (fCurrent==0) fCurrent = fTree->GetFriend(alias);
2992}
2993
2995 TFormLeafInfo(orig)
2996{
2997 // Copy Constructor.
2998 fTree = orig.fTree;
2999 fAlias = orig.fAlias;
3000 fCurrent = orig.fCurrent;
3001}
3002
3004{
3005 // Copy the object and all its content.
3006
3007 return new TFormLeafInfoTTree(*this);
3008}
3009
3010////////////////////////////////////////////////////////////////////////////////
3011/// returns the address of the value pointed to by the
3012/// TFormLeafInfo.
3013
3015{
3016 return GetLocalValuePointer((char*)fCurrent,instance);
3017}
3018
3019////////////////////////////////////////////////////////////////////////////////
3020/// Return result of a leafobject method.
3021
3022template <typename T>
3023T TFormLeafInfoTTree::GetValueImpl(TLeaf *, Int_t instance)
3024{
3025 return ReadTypedValue<T>((char*)fCurrent,instance);
3026}
3027
3028////////////////////////////////////////////////////////////////////////////////
3029/// Return result of a leafobject method.
3030
3031template <typename T>
3032T TFormLeafInfoTTree::ReadValueImpl(char *thisobj, Int_t instance)
3033{
3034 if (fElement) return TFormLeafInfo::ReadTypedValue<T>(thisobj,instance);
3035 else if (fNext) return fNext->ReadTypedValue<T>(thisobj,instance);
3036 else return 0;
3037}
3038
3041
3042////////////////////////////////////////////////////////////////////////////////
3043/// Update after a change of file in a chain
3044
3046{
3047 if (fAlias.Length() && fAlias != fTree->GetName()) {
3049 }
3050 return fCurrent && TFormLeafInfo::Update();
3051}
void Class()
Definition: Class.C:29
ROOT::R::TRInterface & r
Definition: Object.C:4
#define d(i)
Definition: RSha256.hxx:102
unsigned short UShort_t
Definition: RtypesCore.h:36
int Int_t
Definition: RtypesCore.h:41
unsigned char UChar_t
Definition: RtypesCore.h:34
char Char_t
Definition: RtypesCore.h:29
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
unsigned long ULong_t
Definition: RtypesCore.h:51
long Long_t
Definition: RtypesCore.h:50
bool Bool_t
Definition: RtypesCore.h:59
short Short_t
Definition: RtypesCore.h:35
double Double_t
Definition: RtypesCore.h:55
long long Long64_t
Definition: RtypesCore.h:69
unsigned long long ULong64_t
Definition: RtypesCore.h:70
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
EDataType
Definition: TDataType.h:28
@ kNoType_t
Definition: TDataType.h:33
@ kULong64_t
Definition: TDataType.h:32
@ kChar_t
Definition: TDataType.h:29
@ kOther_t
Definition: TDataType.h:32
@ kIsPointer
Definition: TDictionary.h:77
@ kIsReference
Definition: TDictionary.h:81
#define R__ASSERT(e)
Definition: TError.h:96
#define INSTANTIATE_GETVAL(CLASS)
#define INSTANTIATE_READVAL(CLASS)
#define READ_ARRAY(TYPE_t)
#define GET_ARRAY(TYPE_t)
int type
Definition: TGX11.cxx:120
#define gInterpreter
Definition: TInterpreter.h:555
#define gROOT
Definition: TROOT.h:415
char * Form(const char *fmt,...)
#define R__WRITE_LOCKGUARD(mutex)
Array of integers (32 bits per element).
Definition: TArrayI.h:27
void Set(Int_t n)
Set size of this array to n ints.
Definition: TArrayI.cxx:105
Int_t At(Int_t i) const
Definition: TArrayI.h:79
void AddAt(Int_t c, Int_t i)
Add Int_t c at position i. Check for out of bounds.
Definition: TArrayI.cxx:93
void Copy(TArrayI &array) const
Definition: TArrayI.h:42
Int_t GetSize() const
Definition: TArray.h:47
A Branch for the case of an object.
TBranchElement * GetBranchCount() const
Int_t GetID() const
TStreamerInfo * GetInfo() const
Get streamer info for the branch class.
virtual char * GetAddress() const
Get the branch address.
TBranchElement * GetBranchCount2() const
char * GetObject() const
Return a pointer to our object.
Int_t GetNdata() const
A TTree is a list of TBranches.
Definition: TBranch.h:91
TTree * GetTree() const
Definition: TBranch.h:250
TObjArray * GetListOfLeaves()
Definition: TBranch.h:245
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition: TClass.h:75
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:4440
Int_t Size() const
Return size of object of this class.
Definition: TClass.cxx:5454
Bool_t IsLoaded() const
Return true if the shared library of this class is currently in the a process's memory.
Definition: TClass.cxx:5662
Bool_t IsTObject() const
Return kTRUE is the class inherits from TObject.
Definition: TClass.cxx:5688
TVirtualCollectionProxy * GetCollectionProxy() const
Return the proxy describing the collection (if any).
Definition: TClass.cxx:2835
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:2906
An array of clone (identical) objects.
Definition: TClonesArray.h:32
A small helper class to implement casting an object to a different type (equivalent to dynamic_cast)
virtual ~TFormLeafInfoCast()
Destructor.
TString fCastedName
Pointer to the class we are trying to case to.
virtual Int_t GetNdata()
Get the number of element in the entry.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
void Swap(TFormLeafInfoCast &other)
Exception safe swap.
Bool_t fIsTObject
Marked by ReadValue.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
TFormLeafInfoCast(TClass *classptr=0, TClass *casted=0)
Indicated whether the fClass inherits from TObject.
Bool_t fGoodCast
Name of the class we are casting to.
TFormLeafInfoCast & operator=(const TFormLeafInfoCast &orig)
Exception safe assignment operator.
A small helper class to implement reading a data member on a TClonesArray object stored in a TTree.
TFormLeafInfoClones(TClass *classptr=0, Long_t offset=0)
Constructor.
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
virtual Int_t ReadCounterValue(char *where)
Return the current size of the the TClonesArray.
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the current size of the the TClonesArray.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
TFormLeafInfoClones & operator=(const TFormLeafInfoClones &orig)
Exception safe assignement operator.
void Swap(TFormLeafInfoClones &other)
Exception safe swap.
A small helper class to implement reading a data member on a TClonesArray object stored in a TTree.
TFormLeafInfoCollectionObject(TClass *classptr=0, Bool_t fTop=kTRUE)
Constructor.
void Swap(TFormLeafInfoCollectionObject &other)
Exception safe swap.
virtual Double_t ReadValue(char *where, Int_t instance=0)
Return the value of the underlying data member inside the clones array.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the current size of the the TClonesArray.
TFormLeafInfoCollectionObject & operator=(const TFormLeafInfoCollectionObject &orig)
Exception safe assignement operator.
Used to return the size of a collection.
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Not implemented.
TVirtualCollectionProxy * fCollProxy
virtual Double_t ReadValue(char *where, Int_t instance=0)
Return the value of the underlying pointer data member.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all of its content.
TFormLeafInfoCollectionSize()
Constructor.
void Swap(TFormLeafInfoCollectionSize &other)
Exception safe swap.
TFormLeafInfoCollectionSize & operator=(const TFormLeafInfoCollectionSize &orig)
Exception safe assignment operator.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Not implemented.
~TFormLeafInfoCollectionSize()
Destructor.
A small helper class to implement reading a data member on a generic collection object stored in a TT...
TVirtualCollectionProxy * fCollProxy
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
TFormLeafInfoCollection()
Constructor.
virtual TFormLeafInfo * DeepCopy() const
Copy of the object and its content.
virtual Bool_t HasCounter() const
Return true if the underlying data has a array size counter.
void Swap(TFormLeafInfoCollection &other)
Exception safe swap.
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the current size of the the TClonesArray.
virtual Int_t ReadCounterValue(char *where)
Return the size of the underlying array for the current entry in the TTree.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
~TFormLeafInfoCollection()
Destructor.
TStreamerElement * fLocalElement
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
Return the pointer to the clonesArray.
TFormLeafInfoCollection & operator=(const TFormLeafInfoCollection &orig)
Exception safe assignment operator.
A small helper class to implement reading a data member on an object stored in a TTree.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
Return the address of the underlying value.
TFormLeafInfoDirect(TBranchElement *from)
Constructor.
virtual TFormLeafInfo * DeepCopy() const
Copy this object and its content.
virtual Double_t ReadValue(char *, Int_t=0)
Read the value at the given memory location.
Asmall helper class to implement executing a method of an object stored in a TTree.
virtual TClass * GetClass() const
Return the type of the underlying return value.
static TClass * ReturnTClass(TMethodCall *mc)
Return the TClass corresponding to the return type of the function if it is an object type or if the ...
TFormLeafInfoMethod & operator=(const TFormLeafInfoMethod &orig)
Exception safe assignment operator.
~TFormLeafInfoMethod()
Destructor.
virtual Bool_t IsString() const
Return true if the return value is a string.
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 Bool_t IsInteger() const
Return true if the return value is integral.
TMethodCall * fMethod
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
TFormLeafInfoMethod(TClass *classptr=0, TMethodCall *method=0)
Constructor.
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...
A small helper class to implement reading a data member on a variable size array inside a TClonesArra...
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
TFormLeafInfoMultiVarDimClones(TClass *motherclassptr, Long_t offset, TClass *elementclassptr, TFormLeafInfo *parent)
Constructor.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its data.
virtual Double_t GetValue(TLeaf *leaf, Int_t instance=0)
A small helper class to implement reading a data member on a variable size array inside a TClonesArra...
virtual Double_t GetValue(TLeaf *leaf, Int_t instance=0)
TFormLeafInfoMultiVarDimCollection(TClass *motherclassptr, Long_t offset, TClass *elementclassptr, TFormLeafInfo *parent)
Constructor.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
A small helper class to implement reading a data member on a variable size array inside a TClonesArra...
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual Double_t ReadValue(char *, Int_t=0)
Not implemented.
A helper class to implement reading a data member on a variable size array inside a TClonesArray obje...
virtual Int_t GetVarDim()
Return the index of the dimension which varies for each elements of an enclosing array (typically a T...
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
virtual Double_t GetValue(TLeaf *, Int_t=0)
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
virtual Int_t GetSize(Int_t index)
Return the size of the requested sub-array.
virtual void SetPrimaryIndex(Int_t index)
Set the current value of the primary index.
TFormLeafInfoMultiVarDim()
Constructor.
TFormLeafInfo * fCounter2
virtual void SetSecondaryIndex(Int_t index)
Set the current value of the primary index.
void Swap(TFormLeafInfoMultiVarDim &other)
Exception safe swap.
TFormLeafInfoMultiVarDim & operator=(const TFormLeafInfoMultiVarDim &orig)
Exception safe assignment operator.
virtual Int_t GetPrimaryIndex()
Return the index vlaue of the primary index.
virtual void SetSize(Int_t index, Int_t val)
Set the sizes of the sub-array.
virtual Int_t GetSumOfSizes()
Get the total size.
virtual Int_t GetVirtVarDim()
Return the virtual index (for this expression) of the dimension which varies for each elements of an ...
~TFormLeafInfoMultiVarDim()
Destructor.
virtual void UpdateSizes(TArrayI *garr)
Update the sizes of the arrays.
A small helper class to implement reading a numerical value inside a collection.
TFormLeafInfoNumerical & operator=(const TFormLeafInfoNumerical &orig)
Exception safe assignment operator.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
void Swap(TFormLeafInfoNumerical &other)
Exception safe swap.
virtual Bool_t IsString() const
Return true if the underlying data is a string.
TFormLeafInfoNumerical(TVirtualCollectionProxy *holder_of)
Construct a TFormLeafInfo for the numerical type contained in the collection.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all its content.
virtual ~TFormLeafInfoNumerical()
Destructor.
A small helper class to implement reading a data member by following a pointer inside a branch of TTr...
TFormLeafInfoPointer(TClass *classptr=0, Long_t offset=0, TStreamerElement *element=0)
Constructor.
virtual TFormLeafInfo * DeepCopy() const
Copy the object and all of its contnet.
A small helper class to implement reading from the containing TTree object itself.
virtual TFormLeafInfo * DeepCopy() const
Make a complete copy of this FormLeafInfo and all its content.
TFormLeafInfoTTree(TTree *tree=0, const char *alias=0, TTree *current=0)
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
returns the address of the value pointed to by the TFormLeafInfo.
virtual Bool_t Update()
Update after a change of file in a chain.
This class is a small helper class to implement reading a data member on an object stored in a TTree.
Definition: TFormLeafInfo.h:49
virtual Int_t GetVarDim()
Return the index of the dimension which varies for each elements of an enclosing array (typically a T...
virtual void * GetValuePointer(TLeaf *leaf, Int_t instance=0)
returns the address of the value pointed to by the serie of TFormLeafInfo.
TFormLeafInfo * fCounter
Descriptor of the data pointed to.
Definition: TFormLeafInfo.h:68
virtual Int_t GetNdata()
Get the number of element in the entry.
TClass * fClass
Definition: TFormLeafInfo.h:62
virtual ~TFormLeafInfo()
Delete this object and all its content.
TFormLeafInfo(TClass *classptr=0, Long_t offset=0, TStreamerElement *element=0)
Constructor.
char * GetObjectAddress(TLeafElement *leaf, Int_t &instance)
Returns the the location of the object pointed to.
TString fElementName
Definition: TFormLeafInfo.h:71
virtual Int_t GetPrimaryIndex()
Method for multiple variable dimensions.
virtual void SetSize(Int_t index, Int_t val)
Set the current size of the arrays.
virtual Bool_t IsInteger() const
Return true if the underlying data is an integral value.
virtual Bool_t HasCounter() const
Return true if any of underlying data has a array size counter.
virtual TClass * GetClass() const
Get the class of the underlying data.
TStreamerElement * fElement
Offset of the data pointed inside the class fClass.
Definition: TFormLeafInfo.h:65
TFormLeafInfo * fNext
Definition: TFormLeafInfo.h:69
virtual Int_t GetVirtVarDim()
Return the virtual index (for this expression) of the dimension which varies for each elements of an ...
virtual void SetSecondaryIndex(Int_t index)
Set the primary index value.
Long_t fOffset
This is the class of the data pointed to.
Definition: TFormLeafInfo.h:64
virtual void LoadSizes(TBranch *branch)
Load the current array sizes.
virtual Int_t GetSumOfSizes()
Total all the elements that are available for the current entry for the secondary variable dimension.
virtual Int_t GetArrayLength()
Return the current length of the array.
TFormLeafInfo & operator=(const TFormLeafInfo &orig)
Exception safe assignment operator.
T ReadTypedValue(char *where, Int_t instance=0)
void Swap(TFormLeafInfo &other)
virtual Int_t ReadCounterValue(char *where)
Return the size of the underlying array for the current entry in the TTree.
virtual Bool_t Update()
We reloading all cached information in case the underlying class information has changed (for example...
virtual Int_t GetCounterValue(TLeaf *leaf)
Return the size of the underlying array for the current entry in the TTree.
Int_t GetNdata(TLeaf *leaf)
virtual void UpdateSizes(TArrayI *garr)
Set the current sizes of the arrays.
virtual Int_t GetSize(Int_t index)
For the current entry, and the value 'index' for the main array, return the size of the secondary var...
virtual void SetPrimaryIndex(Int_t index)
Set the primary index value.
virtual void AddOffset(Int_t offset, TStreamerElement *element)
Increase the offset of this element.
Int_t GetMultiplicity()
Reminder of the meaning of fMultiplicity:
virtual Bool_t IsString() const
Return true if the underlying data is a string.
virtual void * GetLocalValuePointer(TLeaf *leaf, Int_t instance=0)
returns the address of the value pointed to by the TFormLeafInfo.
TString fClassName
Definition: TFormLeafInfo.h:70
Int_t fMultiplicity
Definition: TFormLeafInfo.h:74
virtual TFormLeafInfo * DeepCopy() const
Make a complete copy of this FormLeafInfo and all its content.
Long_t Property() const
Get property description word. For meaning of bits see EProperty.
Definition: TFunction.cxx:183
const char * GetReturnTypeName() const
Get full type description of function return type, e,g.: "class TDirectory*".
Definition: TFunction.cxx:140
A TLeaf for the general case when using the branches created via a TStreamerInfo (i....
Definition: TLeafElement.h:32
A TLeaf for a general object derived from TObject.
Definition: TLeafObject.h:31
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:49
virtual void * GetValuePointer() const
Definition: TLeaf.h:129
T GetTypedValue(Int_t i=0) const
Definition: TLeaf.h:135
TBranch * GetBranch() const
Definition: TLeaf.h:107
Method or function calling interface.
Definition: TMethodCall.h:37
EReturnType ReturnType()
Returns the return type of the method.
static const EReturnType kLong
Definition: TMethodCall.h:43
const char * GetMethodName() const
Definition: TMethodCall.h:90
static const EReturnType kString
Definition: TMethodCall.h:45
static const EReturnType kOther
Definition: TMethodCall.h:46
const char * GetParams() const
Definition: TMethodCall.h:91
TFunction * GetMethod()
Returns the TMethod describing the method to be executed.
void Init(const TFunction *func)
Initialize the method invocation environment based on the TFunction object.
void Execute(const char *, const char *, int *=0)
Execute method on this object with the given parameter string, e.g.
Definition: TMethodCall.h:64
static const EReturnType kDouble
Definition: TMethodCall.h:44
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:522
TObject * UncheckedAt(Int_t i) const
Definition: TObjArray.h:90
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
const char * GetCountName() const
Int_t GetNewType() const
virtual TClass * GetClassPointer() const
Returns a pointer to the TClass of this element.
Int_t GetArrayLength() const
TClass * GetClass() const
Describe Streamer information for one class version.
Definition: TStreamerInfo.h:43
@ kUChar
Equal to TDataType's kchar.
Int_t GetNewType(Int_t id) const
Basic string class.
Definition: TString.h:131
Ssiz_t Length() const
Definition: TString.h:405
const char * Data() const
Definition: TString.h:364
TString & Append(const char *cs)
Definition: TString.h:559
A TTree represents a columnar dataset.
Definition: TTree.h:72
virtual TTree * GetFriend(const char *) const
Return a pointer to the TTree friend whose name or alias is 'friendname.
Definition: TTree.cxx:5835
Int_t GetMakeClass() const
Definition: TTree.h:482
virtual EDataType GetType() const =0
virtual void * At(UInt_t idx)=0
virtual UInt_t Size() const =0
virtual TVirtualCollectionProxy * Generate() const =0
virtual Bool_t HasPointers() const =0
virtual TClass * GetCollectionClass() const
std::string GetName(const std::string &scope_name)
Definition: Cppyy.cxx:150
void swap(RDirectoryEntry &e1, RDirectoryEntry &e2) noexcept
double T(double x)
Definition: ChebyshevPol.h:34
R__EXTERN TVirtualRWMutex * gCoreMutex
TClass * GetClass(T *)
Definition: TClass.h:608
static Roo_reg_AGKInteg1D instance
void GetNormalizedName(std::string &norm_name, std::string_view name)
Return the normalized name.
Definition: TClassEdit.cxx:832
Definition: tree.py:1
auto * l
Definition: textangle.C:4