1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
19 // Read the old ALICE event format based on CW-ntuples
20 // http://consult.cern.ch/alice/Internal_Notes/1995/32/abstract
21 // .cwn file have to be converted to .root using h2root
22 // Use SetFileName(file) to read from "file"
23 // Author: andreas.morsch@cern.ch
26 #include <TParticle.h>
28 #include <TVirtualMC.h>
30 #include "AliGenReaderCwn.h"
32 ClassImp(AliGenReaderCwn);
35 AliGenReaderCwn::AliGenReaderCwn()
37 // Default constructor
42 AliGenReaderCwn::~AliGenReaderCwn()
47 void AliGenReaderCwn::Init()
50 // reset the existing file environment and open a new root file if
51 // the pointer to the Fluka tree is null
55 pFile = new TFile(fFileName);
57 printf("\n I have opened %s file \n", fFileName);
59 // get the tree address in the Fluka boundary source file
60 fTreeNtuple = (TTree*)gDirectory->Get("h888");
62 TTree *h2=fTreeNtuple;
63 //Set branch addresses
64 h2->SetBranchAddress("Nihead",&fNihead);
65 h2->SetBranchAddress("Ihead",fIhead);
66 h2->SetBranchAddress("Nrhead",&fNrhead);
67 h2->SetBranchAddress("Rhead",fRhead);
68 h2->SetBranchAddress("Idpart",&fIdpart);
69 h2->SetBranchAddress("Theta",&fTheta);
70 h2->SetBranchAddress("Phi",&fPhi);
71 h2->SetBranchAddress("P",&fP);
72 h2->SetBranchAddress("E",&fE);
75 Int_t AliGenReaderCwn::NextEvent()
77 // Read the next event
80 TFile* pFile = fTreeNtuple->GetCurrentFile();
83 Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
84 if (fNcurrent < nentries) {
90 printf("\n This should never happen !\n");
93 printf("\n Next event contains %d tracks! \n", i6);
96 fNparticleMax = nTracks;
103 TParticle* AliGenReaderCwn::NextParticle()
109 // Read the next particle
110 if (fCode == kGEANT3) fIdpart=gMC->PDGFromId(fIdpart);
111 Double_t amass = TDatabasePDG::Instance()->GetParticle(fIdpart)->Mass();
113 Warning("Generate","Particle %d E = %f mass = %f %f %f \n",
114 fIdpart,fE,amass, fPhi, fTheta);
117 prwn=sqrt((fE+amass)*(fE-amass));
120 fTheta *= TMath::Pi()/180.;
121 fPhi = (fPhi-180)*TMath::Pi()/180.;
122 p[0] = prwn*TMath::Sin(fTheta)*TMath::Cos(fPhi);
123 p[1] = prwn*TMath::Sin(fTheta)*TMath::Sin(fPhi);
124 p[2] = prwn*TMath::Cos(fTheta);
126 TParticle* particle = new TParticle(fIdpart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 0., 0., 0., 0.);
134 AliGenReaderCwn& AliGenReaderCwn::operator=(const AliGenReaderCwn& rhs)
136 // Assignment operator
141 void AliGenReaderCwn::Copy(TObject&) const
146 Fatal("Copy","Not implemented!\n");