Re: [ROOT] Tree and TClonesArray Problem

From: Rene Brun (Rene.Brun@cern.ch)
Date: Wed Dec 13 2000 - 14:28:22 MET


Hi Claus,
Change the line
      m_Tree->SetBranchAddress("particles",m_Fruits);
to
      m_Tree->SetBranchAddress("particles",&m_Fruits);

Rene Brun


Claus Peter Buszello wrote:
> 
> Hello,
> 
> I am stuck with a problem trying to fill a tree containing TClonesArrays,
> and reading it again.
> 
> This is how I write the file:
> 
> #include "TROOT.h"
> #include <TFile.h>
> #include <TTree.h>
> #include <TClonesArray.h>
> #include <TMCParticle.h>
> 
> main(){
>   char buf[65535];
>   TROOT api("Test","Test");
> 
>   TFile *F = new TFile("pythia5.root","RECREATE");
>   TTree *tree = new TTree("t","t");
> 
>   TClonesArray * particles = new TClonesArray("TMCParticle",10000,kFALSE);
>   TClonesArray &p = *particles;
>   tree->Branch("particles",&particles,32000,0);
>   while (fgets(buf,65534,stdin) != NULL){
> 
> // unimportant
>   long double P1[4],P2[4],P3[4],P4[4];
>   sscanf(buf,"%*llE %*llE %*llE %*llE %llE %llE %llE %llE %*llE %llE %llE
> %llE %llE %*llE %llE %llE %llE %llE %*llE %llE %llE %llE %llE",
>            &P1[0],&P1[1],&P1[2],&P1[3],
>            &P2[0],&P2[1],&P2[2],&P2[3],
>            &P3[0],&P3[1],&P3[2],&P3[3],
>            &P4[0],&P4[1],&P4[2],&P4[3]);
> 
>     new(p[0]) TMCParticle(1,13,0,0,0,
>                       P1[0],P1[1],P1[2],P1[3],
>                       0.105658,0,0,0,0,0);
>     new(p[1]) TMCParticle(1,-13,0,0,0,
>                       P2[0],P2[1],P2[2],P2[3],
>                       0.105658,0,0,0,0,0);
>     new(p[2]) TMCParticle(1,13,0,0,0,
>                       P3[0],P3[1],P3[2],P3[3],
>                       0.105658,0,0,0,0,0);
>     new(p[3]) TMCParticle(1,-13,0,0,0,
>                       P4[0],P4[1],P4[2],P4[3],
>                       0.105658,0,0,0,0,0);
> // end unimportant
> 
>    tree->Fill();
>   }
>   tree->Print();
>   F->Write();
>   F->Close();
> }
> 
> After this has run I can browse the file and everything looks and feels
> OK.
> 
> Then the tree is read by a code like the following:
> 
> #include "TROOT.h"
> #include <TFile.h>
> #include <TTree.h>
> #include <TClonesArray.h>
> #include <TMCParticle.h>
> 
> main(){
> 
>      TROOT api("Test","Test");
> 
>      TClonesArray *m_Fruits = new TClonesArray("TMCParticle",10000,
> kFALSE);
> 
>       TFile * m_File = new TFile("pythia5.root");
>       TTree* m_Tree;
>       m_Tree = (TTree*)m_File->Get("t");
> 
>       m_Tree->SetBranchAddress("particles",m_Fruits);
>       for (int i=0;i<100;i++){
>          m_Tree->GetEvent(i);
>          printf("%i\n",m_Fruits->GetEntriesFast());
>       }
> }
> 
> The Problem is, that m_Fruits->GetEntriesFast() always returns 0 even
> though there is data in the file.
> 
> The Code that reads the file is the one that is most likely to be allright
> (it is part of atlfast++) and shouldn't be changed. I have to change my
> code so that it fits, and produces Files, that can be read by atlfast++.
> 
> If You have any Idea what I have done wrong I would be thankful if you
> let me know (I have tried a lot and have no more Ideas...)
> 
> Greetings
> 
>    Claus Peter Buszello



This archive was generated by hypermail 2b29 : Tue Jan 02 2001 - 11:50:39 MET