Re: [ROOT] hepevt to root ntuple conversion

From: Rene Brun (Rene.Brun@cern.ch)
Date: Wed Jan 03 2001 - 09:40:21 MET


Hi,
A toy example using hepevt is following. You can also have a look
at TGenerator::ImportParticles to see how to convert a hepevt structure
into a TClonesArray of TParticle objects.

Rene Brun

typedef struct {
	Int_t	 nevhep;
        Int_t    nhep;
        Int_t    isthep[2000];
        Int_t    idhep[2000];
        Int_t    jmohep[2000][2];
        Int_t    jdahep[2000][2];
        Double_t phep[2000][5];
        Double_t vhep[2000][4];
} HEPEVT;

void hepevt()
{
   HEPEVT event;

   TFile f("hepevt.root","recreate");
   TTree *T = new TTree("T","hepevt ntuple");
   T->Branch("nhep",  &event.nhep,  "nhep/I");
   T->Branch("isthep",&event.isthep,"isthep[nhep]/I");
   T->Branch("idhep", &event.idhep, "idhep[nhep]/I");
   T->Branch("jmohep",&event.jmohep,"jmohep[nhep][2]/I");
   T->Branch("jdahep",&event.jdahep,"jdahep[nhep][2]/I");
   T->Branch("phep",  &event.phep,  "phep[nhep][5]/I");
   T->Branch("vhep",  &event.vhep,  "vhep[nhep][4]/I");
   
   TRandom r;
   Double_t px,py,pz,p,E;
   Int_t id, nevent = 100;
   for (Int_t i=0;i<nevent;i++) {
      event.nhep = gRandom->Rndm()*499;
      for (Int_t n=0;n<event.nhep;n++) {
         event.isthep[n] = (Int_t)(r.Rndm()*2);
         event.idhep[n]  = id = (Int_t)(r.Rndm()*100);
         event.jmohep[n][0]  = 1;
         event.jmohep[n][1]  = 2;
         event.jdahep[n][0]  = 0;
         event.jdahep[n][1]  = 0;
         event.phep[n][0]  = px = r.Gaus(0,1);
         event.phep[n][1]  = py = r.Gaus(0,1);
         event.phep[n][2]  = pz = r.Gaus(0,10);
         event.phep[n][3]  = p  = TMath::Sqrt(px*px+py*py+pz*pz);
         event.vhep[n][4]  = E  = TMath::Sqrt(p*p+id);
         event.vhep[n][0]  = r.Gaus(0,1e-3);
         event.vhep[n][1]  = r.Gaus(0,1e-3);
         event.vhep[n][2]  = r.Gaus(0,10);
         event.vhep[n][3]  = r.Gaus(0,1e-9);
      }
      T->Fill();
   }
   T->Print();
   T->Write();
   f.Close();
}


On Tue, 2 Jan 2001, Muge Karagoz wrote:

> Hello all,
> 
> I have a file in hepevt format.I would like to convert the contents of
> this file to a root ntuple.
> Is there an example/already existing program to do it? Could anyone please
> give me a suggestion?
> 
> Thank you very much,
> 
> 
> 



This archive was generated by hypermail 2b29 : Tue Jan 01 2002 - 17:50:33 MET