Re: [ROOT] variable length branches in TTree's

From: Rene Brun (Rene.Brun@cern.ch)
Date: Fri Jun 16 2000 - 11:24:32 MEST


Hi Peter,
If you do not want to use a class (see Pasha's answer), you may get some
inspiration from this simple example:

Rene Brun

{
   TFile *f = new TFile("peter.root","recreate");
   Int_t nPhot;
   Float_t E[500];
   Float_t Px[500];
   Float_t Py[500];
   Float_t Pz[500];
   Float_t Arm[500];
   Float_t PhotProb[500];

   TTree* nEmcPhotons = new TTree("nEmcPhotons","EMC Photons");
   nEmcPhotons->Branch("nPhot",&nPhot,"nPhot/I");
   nEmcPhotons->Branch("E",E,"E[nPhot]/F");
   nEmcPhotons->Branch("Px",Px,"Px[nPhot]/F");
   nEmcPhotons->Branch("Py",Py,"Py[nPhot]/F");
   nEmcPhotons->Branch("Pz",Pz,"Pz[nPhot]/F");
   nEmcPhotons->Branch("Arm", Arm,"Arm[nPhot]/F");
   nEmcPhotons->Branch("PhotProb", PhotProb,"PhotProb[nPhot]/F");

   for (Int_t i=0;i<1000;i++) {
      nPhot = 500*gRandom->Rndm();
      for (Int_t j=0;j<nPhot;j++) {
         E[j] = j;
         Px[j] = j+1;
         Py[j] = j+2;
         Pz[j] = j+3;
         Arm[j] = j+4;
         PhotProb[j] = j+5;
      }
      nEmcPhotons->Fill();
   }
   nEmcPhotons->Print();
   nEmcPhotons->Write();
   f->Close();
   delete f;
}

Peter Steinberg wrote:
> 
> Hi -
> 
> Is it straightforward to have variable-length branches
> in a TTree?  I would like to have functionality a la HBOOK
> CWN's, i.e. a branch which has a leaflist like:
> 
> "N:x[N]/F:y[N]/F:z[N]/F" etc.
> 
> where i point the branch to a structure defined somewhere:
> 
> {
>   Int_t N;
>   Float_t x[32];
>   Float_t y[32];
>   Float_t z[32];
> }
> 
> The point is that when I set N to 10, then only 10 x,y,z's will be
> filled in.
> 
> I am somewhat of a novice with trees, so this may be simple, but
> i have been unable to make anything work.
> 
> Regards,
> Peter
> 
> ---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+
> Dr. Peter A. Steinberg                        Brookhaven National Laboratory
> Assistant Scientist                                Bldg 555, Upton NY, 11973
> ---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+
> Office: 631-344-4390                          mailto:Peter.Steinberg@bnl.gov
> Mobile: 917-549-3094                       http://www.rhic.bnl.gov/~steinber
> --> E-mail my cell phone (100 chars max): steinber@steinberg.chm.bnl.gov <--
> ---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+---+



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