Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TFormLeafInfo.h
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 #ifndef ROOT_TFormLeafInfo
13 #define ROOT_TFormLeafInfo
14 
15 #include "TObject.h"
16 
17 #include "TLeafElement.h"
18 
19 #include "TArrayI.h"
20 #include "TDataType.h"
21 #include "TStreamerInfo.h"
22 #include "TStreamerElement.h"
23 
24 
25 // declare the extra versions of GetValue() plus templated implementation
26 #define DECLARE_GETVAL \
27  virtual Double_t GetValue(TLeaf *leaf, Int_t instance = 0) \
28  { return GetValueImpl<Double_t>(leaf, instance); } \
29  virtual Long64_t GetValueLong64(TLeaf *leaf, Int_t instance = 0) \
30  { return GetValueImpl<Long64_t>(leaf, instance); } \
31  virtual LongDouble_t GetValueLongDouble(TLeaf *leaf, Int_t instance = 0) \
32  { return GetValueImpl<LongDouble_t>(leaf, instance); } \
33  template<typename T> T GetValueImpl(TLeaf *leaf, Int_t instance = 0) // no semicolon
34 
35 
36 // declare the extra versions of ReadValue() plus templated implementation
37 #define DECLARE_READVAL \
38  virtual Double_t ReadValue(char *where, Int_t instance = 0) \
39  { return ReadValueImpl<Double_t>(where, instance); } \
40  virtual Long64_t ReadValueLong64(char *where, Int_t instance = 0) \
41  { return ReadValueImpl<Long64_t>(where, instance); } \
42  virtual LongDouble_t ReadValueLongDouble(char *where, Int_t instance = 0) \
43  { return ReadValueImpl<LongDouble_t>(where, instance); } \
44  template<typename T> T ReadValueImpl(char *where, Int_t instance = 0) // no semicolon
45 
46 
47 
48 
49 class TFormLeafInfo : public TObject {
50 public:
51  // Constructors
52  TFormLeafInfo(TClass* classptr = 0, Long_t offset = 0,
53  TStreamerElement* element = 0);
54  TFormLeafInfo(const TFormLeafInfo& orig);
55  virtual TFormLeafInfo* DeepCopy() const;
56  virtual ~TFormLeafInfo();
57 
58  void Swap(TFormLeafInfo &other);
59  TFormLeafInfo &operator=(const TFormLeafInfo &orig);
60 
61  // Data Members
62  TClass *fClass; //! This is the class of the data pointed to
63  //TStreamerInfo *fInfo; //! == fClass->GetStreamerInfo()
64  Long_t fOffset; //! Offset of the data pointed inside the class fClass
65  TStreamerElement *fElement; //! Descriptor of the data pointed to.
66  //Warning, the offset in fElement is NOT correct because it does not take into
67  //account base classes and nested objects (which fOffset does).
68  TFormLeafInfo *fCounter;
69  TFormLeafInfo *fNext; // follow this to grab the inside information
70  TString fClassName;
71  TString fElementName;
72 
73 protected:
74  Int_t fMultiplicity;
75 public:
76 
77  virtual void AddOffset(Int_t offset, TStreamerElement* element);
78 
79  virtual Int_t GetArrayLength();
80  virtual TClass* GetClass() const;
81  virtual Int_t GetCounterValue(TLeaf* leaf);
82  virtual Int_t ReadCounterValue(char *where);
83 
84  char* GetObjectAddress(TLeafElement* leaf, Int_t &instance);
85 
86  Int_t GetMultiplicity();
87 
88  // Currently only implemented in TFormLeafInfoCast
89  Int_t GetNdata(TLeaf* leaf);
90  virtual Int_t GetNdata();
91 
92  virtual void *GetValuePointer(TLeaf *leaf, Int_t instance = 0);
93  virtual void *GetValuePointer(char *from, Int_t instance = 0);
94  virtual void *GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0);
95  virtual void *GetLocalValuePointer( char *from, Int_t instance = 0);
96 
97  virtual Bool_t HasCounter() const;
98  virtual Bool_t IsString() const;
99 
100  virtual Bool_t IsInteger() const;
101  virtual Bool_t IsReference() const { return kFALSE; }
102 
103  // Method for multiple variable dimensions.
104  virtual Int_t GetPrimaryIndex();
105  virtual Int_t GetVarDim();
106  virtual Int_t GetVirtVarDim();
107  virtual Int_t GetSize(Int_t index);
108  virtual Int_t GetSumOfSizes();
109  virtual void LoadSizes(TBranch* branch);
110  virtual void SetPrimaryIndex(Int_t index);
111  virtual void SetSecondaryIndex(Int_t index);
112  virtual void SetSize(Int_t index, Int_t val);
113  virtual void SetBranch(TBranch* br) { if ( fNext ) fNext->SetBranch(br); }
114  virtual void UpdateSizes(TArrayI *garr);
115 
116  virtual Bool_t Update();
117 
118  DECLARE_GETVAL;
119  DECLARE_READVAL;
120 
121  template <typename T> struct ReadValueHelper {
122  static T Exec(TFormLeafInfo *leaf, char *where, Int_t instance) {
123  return leaf->ReadValue(where, instance);
124  }
125  };
126  template <typename T > T ReadTypedValue(char *where, Int_t instance = 0) {
127  return ReadValueHelper<T>::Exec(this, where, instance);
128  }
129 
130  template <typename T> struct GetValueHelper {
131  static T Exec(TFormLeafInfo *linfo, TLeaf *leaf, Int_t instance) {
132  return linfo->GetValue(leaf, instance);
133  }
134  };
135  template <typename T > T GetTypedValue(TLeaf *leaf, Int_t instance = 0) {
136  return GetValueHelper<T>::Exec(this, leaf, instance);
137  }
138 };
139 
140 
141 template <> struct TFormLeafInfo::ReadValueHelper<Long64_t> {
142  static Long64_t Exec(TFormLeafInfo *leaf, char *where, Int_t instance) { return leaf->ReadValueLong64(where, instance); }
143 };
144 template <> struct TFormLeafInfo::ReadValueHelper<ULong64_t> {
145  static ULong64_t Exec(TFormLeafInfo *leaf, char *where, Int_t instance) { return (ULong64_t)leaf->ReadValueLong64(where, instance); }
146 };
147 template <> struct TFormLeafInfo::ReadValueHelper<LongDouble_t> {
148  static LongDouble_t Exec(TFormLeafInfo *leaf, char *where, Int_t instance) { return leaf->ReadValueLongDouble(where, instance); }
149 };
150 
151 template <> struct TFormLeafInfo::GetValueHelper<Long64_t> {
152  static Long64_t Exec(TFormLeafInfo *linfo, TLeaf *leaf, Int_t instance) { return linfo->GetValueLong64(leaf, instance); }
153 };
154 template <> struct TFormLeafInfo::GetValueHelper<ULong64_t> {
155  static ULong64_t Exec(TFormLeafInfo *linfo, TLeaf *leaf, Int_t instance) { return (ULong64_t)linfo->GetValueLong64(leaf, instance); }
156 };
157 template <> struct TFormLeafInfo::GetValueHelper<LongDouble_t> {
158  static LongDouble_t Exec(TFormLeafInfo *linfo, TLeaf *leaf, Int_t instance) { return linfo->GetValueLongDouble(leaf, instance); }
159 };
160 
161 // TFormLeafInfoDirect is a small helper class to implement reading a data
162 // member on an object stored in a TTree.
163 
164 class TFormLeafInfoDirect : public TFormLeafInfo {
165 public:
166  TFormLeafInfoDirect(TBranchElement * from);
167  // The implicit default constructor's implementation is correct.
168 
169  virtual TFormLeafInfo* DeepCopy() const;
170 
171  DECLARE_GETVAL;
172  virtual void *GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0);
173  virtual void *GetLocalValuePointer(char *thisobj, Int_t instance = 0);
174 
175  virtual Double_t ReadValue(char * /*where*/, Int_t /*instance*/= 0);
176  virtual Long64_t ReadValueLong64(char *where, Int_t i= 0) { return ReadValue(where, i); }
177  virtual LongDouble_t ReadValueLongDouble(char *where, Int_t i= 0) { return ReadValue(where, i); }
178 
179 };
180 
181 // TFormLeafInfoNumerical is a small helper class to implement reading a
182 // numerical value inside a collection
183 
184 class TFormLeafInfoNumerical : public TFormLeafInfo {
185  EDataType fKind;
186  Bool_t fIsBool;
187 public:
188  TFormLeafInfoNumerical(TVirtualCollectionProxy *holder_of);
189  TFormLeafInfoNumerical(EDataType kind);
190  TFormLeafInfoNumerical(const TFormLeafInfoNumerical& orig);
191 
192  virtual TFormLeafInfo* DeepCopy() const;
193  void Swap(TFormLeafInfoNumerical &other);
194  TFormLeafInfoNumerical &operator=(const TFormLeafInfoNumerical &orig);
195 
196  virtual ~TFormLeafInfoNumerical();
197 
198  virtual Bool_t IsString() const;
199  virtual Bool_t Update();
200 };
201 
202 // TFormLeafInfoCollectionObject
203 // This class is used when we are interested by the collection it self and
204 // it is split.
205 
206 class TFormLeafInfoCollectionObject : public TFormLeafInfo {
207  Bool_t fTop; //If true, it indicates that the branch itself contains
208 public:
209  TFormLeafInfoCollectionObject(TClass* classptr = 0, Bool_t fTop = kTRUE);
210  TFormLeafInfoCollectionObject(const TFormLeafInfoCollectionObject &orig);
211 
212  void Swap(TFormLeafInfoCollectionObject &other);
213  TFormLeafInfoCollectionObject &operator=(const TFormLeafInfoCollectionObject &orig);
214 
215  virtual TFormLeafInfo* DeepCopy() const {
216  return new TFormLeafInfoCollectionObject(*this);
217  }
218 
219  DECLARE_GETVAL;
220  virtual Int_t GetCounterValue(TLeaf* leaf);
221  virtual Double_t ReadValue(char *where, Int_t instance = 0);
222  virtual Long64_t ReadValueLong64(char *where, Int_t i= 0) { return ReadValue(where, i); }
223  virtual LongDouble_t ReadValueLongDouble(char *where, Int_t i= 0) { return ReadValue(where, i); }
224  virtual void *GetValuePointer(TLeaf *leaf, Int_t instance = 0);
225  virtual void *GetValuePointer(char *thisobj, Int_t instance = 0);
226  virtual void *GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0);
227  virtual void *GetLocalValuePointer(char *thisobj, Int_t instance = 0);
228 };
229 
230 // TFormLeafInfoClones is a small helper class to implement reading a data
231 // member on a TClonesArray object stored in a TTree.
232 
233 class TFormLeafInfoClones : public TFormLeafInfo {
234  Bool_t fTop; //If true, it indicates that the branch itself contains
235 public:
236  //either the clonesArrays or something inside the clonesArray
237  TFormLeafInfoClones(TClass* classptr = 0, Long_t offset = 0);
238  TFormLeafInfoClones(TClass* classptr, Long_t offset, Bool_t top);
239  TFormLeafInfoClones(TClass* classptr, Long_t offset, TStreamerElement* element,
240  Bool_t top = kFALSE);
241  TFormLeafInfoClones(const TFormLeafInfoClones &orig);
242 
243  void Swap(TFormLeafInfoClones &other);
244  TFormLeafInfoClones &operator=(const TFormLeafInfoClones &orig);
245 
246  virtual TFormLeafInfo* DeepCopy() const {
247  return new TFormLeafInfoClones(*this);
248  }
249 
250  DECLARE_GETVAL;
251  DECLARE_READVAL;
252  virtual Int_t GetCounterValue(TLeaf* leaf);
253  virtual Int_t ReadCounterValue(char *where);
254  virtual void *GetValuePointer(TLeaf *leaf, Int_t instance = 0);
255  virtual void *GetValuePointer(char *thisobj, Int_t instance = 0);
256  virtual void *GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0);
257  virtual void *GetLocalValuePointer(char *thisobj, Int_t instance = 0);
258 };
259 
260 // TFormLeafInfoCollection is a small helper class to implement reading a data member
261 // on a generic collection object stored in a TTree.
262 
263 class TFormLeafInfoCollection : public TFormLeafInfo {
264  Bool_t fTop; //If true, it indicates that the branch itself contains
265  //either the clonesArrays or something inside the clonesArray
266  TClass *fCollClass;
267  TString fCollClassName;
268  TVirtualCollectionProxy *fCollProxy;
269  TStreamerElement *fLocalElement;
270 public:
271 
272  TFormLeafInfoCollection(TClass* classptr,
273  Long_t offset,
274  TStreamerElement* element,
275  Bool_t top = kFALSE);
276 
277  TFormLeafInfoCollection(TClass* motherclassptr,
278  Long_t offset = 0,
279  TClass* elementclassptr = 0,
280  Bool_t top = kFALSE);
281 
282  TFormLeafInfoCollection();
283  TFormLeafInfoCollection(const TFormLeafInfoCollection& orig);
284 
285  ~TFormLeafInfoCollection();
286 
287  void Swap(TFormLeafInfoCollection &other);
288  TFormLeafInfoCollection &operator=(const TFormLeafInfoCollection &orig);
289 
290  virtual TFormLeafInfo* DeepCopy() const;
291 
292  virtual Bool_t Update();
293 
294  DECLARE_GETVAL;
295  DECLARE_READVAL;
296  virtual Int_t GetCounterValue(TLeaf* leaf);
297  virtual Int_t ReadCounterValue(char* where);
298  virtual Int_t GetCounterValue(TLeaf* leaf, Int_t instance);
299  virtual Bool_t HasCounter() const;
300  virtual void *GetValuePointer(TLeaf *leaf, Int_t instance = 0);
301  virtual void *GetValuePointer(char *thisobj, Int_t instance = 0);
302  virtual void *GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0);
303  virtual void *GetLocalValuePointer(char *thisobj, Int_t instance = 0);
304 };
305 
306 // TFormLeafInfoCollectionSize is used to return the size of a collection
307 
308 class TFormLeafInfoCollectionSize : public TFormLeafInfo {
309  TClass *fCollClass;
310  TString fCollClassName;
311  TVirtualCollectionProxy *fCollProxy;
312 public:
313  TFormLeafInfoCollectionSize(TClass*);
314  TFormLeafInfoCollectionSize(TClass* classptr,Long_t offset,TStreamerElement* element);
315  TFormLeafInfoCollectionSize();
316  TFormLeafInfoCollectionSize(const TFormLeafInfoCollectionSize& orig);
317 
318  ~TFormLeafInfoCollectionSize();
319 
320  void Swap(TFormLeafInfoCollectionSize &other);
321  TFormLeafInfoCollectionSize &operator=(const TFormLeafInfoCollectionSize &orig);
322 
323  virtual TFormLeafInfo* DeepCopy() const;
324 
325  virtual Bool_t Update();
326 
327  virtual void *GetValuePointer(TLeaf *leaf, Int_t instance = 0);
328  virtual void *GetValuePointer(char *from, Int_t instance = 0);
329  virtual void *GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0);
330  virtual void *GetLocalValuePointer( char *from, Int_t instance = 0);
331  virtual Double_t ReadValue(char *where, Int_t instance = 0);
332  virtual Long64_t ReadValueLong64(char *where, Int_t i= 0) { return ReadValue(where, i); }
333  virtual LongDouble_t ReadValueLongDouble(char *where, Int_t i= 0) { return ReadValue(where, i); }
334 };
335 
336 // TFormLeafInfoPointer is a small helper class to implement reading a data
337 // member by following a pointer inside a branch of TTree.
338 
339 class TFormLeafInfoPointer : public TFormLeafInfo {
340 public:
341  TFormLeafInfoPointer(TClass* classptr = 0, Long_t offset = 0,
342  TStreamerElement* element = 0);
343  // The default copy constructor is the right implementation.
344 
345  virtual TFormLeafInfo* DeepCopy() const;
346 
347  DECLARE_GETVAL;
348  DECLARE_READVAL;
349 };
350 
351 // TFormLeafInfoMethod is a small helper class to implement executing a method
352 // of an object stored in a TTree
353 
354 class TFormLeafInfoMethod : public TFormLeafInfo {
355 
356  TMethodCall *fMethod;
357  TString fMethodName;
358  TString fParams;
359  Double_t fResult;
360  TString fCopyFormat;
361  TString fDeleteFormat;
362  void *fValuePointer;
363  Bool_t fIsByValue;
364 
365 public:
366  static TClass *ReturnTClass(TMethodCall *mc);
367 
368  TFormLeafInfoMethod(TClass* classptr = 0, TMethodCall *method = 0);
369  TFormLeafInfoMethod(const TFormLeafInfoMethod& orig);
370  ~TFormLeafInfoMethod();
371 
372  void Swap(TFormLeafInfoMethod &other);
373  TFormLeafInfoMethod &operator=(const TFormLeafInfoMethod &orig);
374 
375  virtual TFormLeafInfo* DeepCopy() const;
376 
377  DECLARE_READVAL;
378  virtual TClass* GetClass() const;
379  virtual void *GetLocalValuePointer( TLeaf *from, Int_t instance = 0);
380  virtual void *GetLocalValuePointer(char *from, Int_t instance = 0);
381  virtual Bool_t IsInteger() const;
382  virtual Bool_t IsString() const;
383  virtual Bool_t Update();
384 };
385 
386 // TFormLeafInfoMultiVarDim is a small helper class to implement reading a
387 // data member on a variable size array inside a TClonesArray object stored in
388 // a TTree. This is the version used when the data member is inside a
389 // non-split object.
390 
391 class TFormLeafInfoMultiVarDim : public TFormLeafInfo {
392 public:
393  Int_t fNsize;
394  TArrayI fSizes; // Array of sizes of the variable dimension
395  TFormLeafInfo *fCounter2; // Information on how to read the secondary dimensions
396  Int_t fSumOfSizes; // Sum of the content of fSizes
397  Int_t fDim; // physical number of the dimension that is variable
398  Int_t fVirtDim; // number of the virtual dimension to which this object correspond.
399  Int_t fPrimaryIndex; // Index of the dimensions that is indexing the second dimension's size
400  Int_t fSecondaryIndex; // Index of the second dimension
401 
402 protected:
403  TFormLeafInfoMultiVarDim(TClass* classptr, Long_t offset,
404  TStreamerElement* element) : TFormLeafInfo(classptr,offset,element),fNsize(0),fSizes(),fCounter2(0),fSumOfSizes(0),fDim(0),fVirtDim(0),fPrimaryIndex(-1),fSecondaryIndex(-1) {}
405 
406 public:
407  TFormLeafInfoMultiVarDim(TClass* classptr, Long_t offset,
408  TStreamerElement* element, TFormLeafInfo* parent);
409  TFormLeafInfoMultiVarDim();
410  TFormLeafInfoMultiVarDim(const TFormLeafInfoMultiVarDim& orig);
411  ~TFormLeafInfoMultiVarDim();
412 
413  void Swap(TFormLeafInfoMultiVarDim &other);
414  TFormLeafInfoMultiVarDim &operator=(const TFormLeafInfoMultiVarDim &orig);
415 
416  virtual TFormLeafInfo* DeepCopy() const;
417 
418  /* The proper indexing and unwinding of index is done by prior leafinfo in the chain. */
419  //virtual Double_t ReadValue(char *where, Int_t instance = 0) {
420  // return TFormLeafInfo::ReadValue(where,instance);
421  //}
422 
423  virtual void LoadSizes(TBranch* branch);
424  virtual Int_t GetPrimaryIndex();
425  virtual void SetPrimaryIndex(Int_t index);
426  virtual void SetSecondaryIndex(Int_t index);
427  virtual void SetSize(Int_t index, Int_t val);
428  virtual Int_t GetSize(Int_t index);
429  virtual Int_t GetSumOfSizes();
430  virtual Double_t GetValue(TLeaf * /*leaf*/, Int_t /*instance*/ = 0);
431  virtual Long64_t GetValueLong64(TLeaf *leaf, Int_t i= 0) { return GetValue(leaf, i); }
432  virtual LongDouble_t GetValueLongDouble(TLeaf *leaf, Int_t i= 0) { return GetValue(leaf, i); }
433  virtual Int_t GetVarDim();
434  virtual Int_t GetVirtVarDim();
435  virtual Bool_t Update();
436  virtual void UpdateSizes(TArrayI *garr);
437 };
438 
439 // TFormLeafInfoMultiVarDimDirect is a small helper class to implement reading
440 // a data member on a variable size array inside a TClonesArray object stored
441 // in a TTree. This is the version used for split access
442 
443 class TFormLeafInfoMultiVarDimDirect : public TFormLeafInfoMultiVarDim {
444 public:
445  // The default constructor are the correct implementation.
446 
447  virtual TFormLeafInfo* DeepCopy() const;
448 
449  DECLARE_GETVAL;
450  virtual Double_t ReadValue(char * /*where*/, Int_t /*instance*/ = 0);
451  virtual Long64_t ReadValueLong64(char *where, Int_t i= 0) { return ReadValue(where, i); }
452  virtual LongDouble_t ReadValueLongDouble(char *where, Int_t i= 0) { return ReadValue(where, i); }
453 };
454 
455 // TFormLeafInfoMultiVarDimCollection is a small helper class to implement reading
456 // a data member which is a collection inside a TClonesArray or collection object
457 // stored in a TTree. This is the version used for split access
458 //
459 class TFormLeafInfoMultiVarDimCollection : public TFormLeafInfoMultiVarDim {
460 public:
461  TFormLeafInfoMultiVarDimCollection(TClass* motherclassptr, Long_t offset,
462  TClass* elementclassptr, TFormLeafInfo *parent);
463  TFormLeafInfoMultiVarDimCollection(TClass* classptr, Long_t offset,
464  TStreamerElement* element, TFormLeafInfo* parent);
465  // The default copy constructor is the right implementation.
466 
467  virtual TFormLeafInfo* DeepCopy() const;
468 
469  virtual Int_t GetArrayLength() { return 0; }
470  virtual void LoadSizes(TBranch* branch);
471  virtual Double_t GetValue(TLeaf *leaf, Int_t instance = 0);
472  virtual Long64_t GetValueLong64(TLeaf *leaf, Int_t i= 0) { return GetValue(leaf, i); }
473  virtual LongDouble_t GetValueLongDouble(TLeaf *leaf, Int_t i= 0) { return GetValue(leaf, i); }
474  DECLARE_READVAL;
475 };
476 
477 // TFormLeafInfoMultiVarDimClones is a small helper class to implement reading
478 // a data member which is a TClonesArray inside a TClonesArray or collection object
479 // stored in a TTree. This is the version used for split access
480 //
481 class TFormLeafInfoMultiVarDimClones : public TFormLeafInfoMultiVarDim {
482 public:
483  TFormLeafInfoMultiVarDimClones(TClass* motherclassptr, Long_t offset,
484  TClass* elementclassptr, TFormLeafInfo *parent);
485  TFormLeafInfoMultiVarDimClones(TClass* classptr, Long_t offset,
486  TStreamerElement* element, TFormLeafInfo* parent);
487  // The default copy constructor is the right implementation.
488 
489  virtual TFormLeafInfo* DeepCopy() const;
490 
491  virtual Int_t GetArrayLength() { return 0; }
492  virtual void LoadSizes(TBranch* branch);
493  virtual Double_t GetValue(TLeaf *leaf, Int_t instance = 0);
494  virtual Long64_t GetValueLong64(TLeaf *leaf, Int_t i= 0) { return GetValue(leaf, i); }
495  virtual LongDouble_t GetValueLongDouble(TLeaf *leaf, Int_t i= 0) { return GetValue(leaf, i); }
496  DECLARE_READVAL;
497 };
498 
499 // TFormLeafInfoCast is a small helper class to implement casting an object to
500 // a different type (equivalent to dynamic_cast)
501 
502 class TFormLeafInfoCast : public TFormLeafInfo {
503 public:
504  TClass *fCasted; //! Pointer to the class we are trying to case to
505  TString fCastedName; //! Name of the class we are casting to.
506  Bool_t fGoodCast; //! Marked by ReadValue.
507  Bool_t fIsTObject; //! Indicated whether the fClass inherits from TObject.
508 
509  TFormLeafInfoCast(TClass* classptr = 0, TClass* casted = 0);
510  TFormLeafInfoCast(const TFormLeafInfoCast& orig);
511  virtual ~TFormLeafInfoCast();
512 
513  void Swap(TFormLeafInfoCast &other);
514  TFormLeafInfoCast &operator=(const TFormLeafInfoCast &orig);
515 
516  virtual TFormLeafInfo* DeepCopy() const;
517 
518  DECLARE_READVAL;
519  // Currently only implemented in TFormLeafInfoCast
520  virtual Int_t GetNdata();
521  virtual Bool_t Update();
522 };
523 
524 // TFormLeafInfoTTree is a small helper class to implement reading
525 // from the containing TTree object itself.
526 
527 class TFormLeafInfoTTree : public TFormLeafInfo {
528  TTree *fTree;
529  TTree *fCurrent;
530  TString fAlias;
531 
532 public:
533  TFormLeafInfoTTree(TTree *tree = 0, const char *alias = 0, TTree *current = 0);
534  TFormLeafInfoTTree(const TFormLeafInfoTTree& orig);
535 
536  void Swap(TFormLeafInfoTTree &other);
537  TFormLeafInfoTTree &operator=(const TFormLeafInfoTTree &orig);
538 
539  virtual TFormLeafInfo* DeepCopy() const;
540 
541  using TFormLeafInfo::GetLocalValuePointer;
542  using TFormLeafInfo::GetValue;
543 
544  DECLARE_GETVAL;
545  DECLARE_READVAL;
546  virtual void *GetLocalValuePointer(TLeaf *leaf, Int_t instance = 0);
547  virtual Bool_t Update();
548 };
549 
550 
551 #endif /* ROOT_TFormLeafInfo */
552