Re: [ROOT] TTree::SetMakeClass()

From: Rene Brun (Rene.Brun@cern.ch)
Date: Mon May 19 2003 - 18:12:09 MEST


Hi Vincent,

When the Tree has been created with a top level branch,
you must set the address of this top level branch when you read the 
Tree. You do not have to set the address of individual branches.
By default TTree::GetEntry will rebuilt in memory your TEvent class.

If you are not interested by the TEvent object, but only by one or a few
branches, you have to indicate this intended behaviour to the tree by 
calling TTree::SetMakeClass. The call to SetMakeClass is automatically
generated when generating the analysis code via TTree::MakeClass.

Note that you cannot mix the two modes. If you use SetMakeClass(1),
then decide to move to the normal analysis mode (via TBrowser for 
example), you must call tree.SetMakeClass(0).

Rene Brun


On Mon, 19 May 
2003, Vincent [iso-8859-15] Régnard wrote:

> -----BEGIN PGP SIGNED MESSAGE-----
> Hash: SHA1
> 
> Hi rooters,
> I have been seekink for a bug in my code for hours and I finally found the 
> solution in the roottalkdigest ;) , but would like to know more about it.
> I have some classes similar to the the exemple TEvent. I built a tree with it 
> and now want to access the branches. Drawing, cutting and using the 
> treeviewer did worked pefectly. I am now using the possibility to only read 
> one branch as specified in the documentation:
> 
> - ------- my tree names input_tree  
> // only get kin energy branch
>   input_tree->SetBranchStatus("*",0);// disbles all branches
>   input_tree->SetBranchStatus("Ek",1);// enables kin energy branche
> 
>   TBranch *  Ek_branch=input_tree->GetBranch("Ek");// or get 
>   Ek_branch->SetBasketSize(16000);
>   Ek_branch->SetAddress(&EvKinEnergy);
>   if (!Ek_branch) cout << "Kinetic energy branch undefined" << endl;
>   if (!input_tree->GetBranchStatus("Ek")) cout << "Kinetic energy branch not 
> active" << endl;
>   cout << "Entries in event branch " << Ek_branch->GetEntries() << endl;
>   cout << "Get entry err code : " << Ek_branch->GetEntry(5) << endl;
>   cout << "Energy = " << EvKinEnergy << endl;
> - -----
> 
> when doing that I always had the same value for EvKinEnergy (its initial 
> value) instead of the value stored in the tree. When calling the 
> TBranch::GetEntries() succeded, but reading the values went wrong. I saw that 
> one is to invoque the TTree::SetMakeClass() method to get it work properly:
>   input_tree->SetMakeClass(1);
> but I do not find information on this method and do not hunderstand what it 
> does exactly. Could you please tell me more about the role of this method.
> 
> Thanks a lot for your help.
> 
> - -- 
> V Régnard 
> GPG signature and ID: http://terrier.homelinux.net/html/persopages/pgp.html
> Key fingerprint = F350 E39E 8C1B 4E8E D078  6D77 81D8 B38A E8F3 ED1E
> Key ID: E8F3ED1E
> La vente forcée de logiciels est illégale en France, comment l'empêcher:
> http://www.linuxfrench.net/oem/communique/communique1.html
> -----BEGIN PGP SIGNATURE-----
> Version: GnuPG v1.2.1 (GNU/Linux)
> 
> iD8DBQE+yP37gdiziujz7R4RAjd0AJ0dyVoOMenM8Zs1lwNxyTKCoBta8ACfQSw0
> g2rhdC9ypuvRsST6LDIct7E=
> =AZDN
> -----END PGP SIGNATURE-----
> 



This archive was generated by hypermail 2b29 : Thu Jan 01 2004 - 17:50:11 MET