]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenReaderCwn.cxx
Remove compilation of grndmq
[u/mrichter/AliRoot.git] / EVGEN / AliGenReaderCwn.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16
17 /*
18 $Log$
19 Revision 1.1  2001/11/09 09:10:46  morsch
20 Realisation of AliGenReader that reads the old cwn event format.
21
22 */
23
24 // Read the old ALICE event format based on CW-ntuples
25 // http://consult.cern.ch/alice/Internal_Notes/1995/32/abstract
26 // .cwn file have to be converted to .root using h2root
27 // Use SetFileName(file) to read from "file" 
28 // Author: andreas.morsch@cern.ch
29
30 #include <TFile.h>
31 #include <TTree.h>
32 #include <TParticle.h>
33
34 #include "AliGenReaderCwn.h"
35 #include "AliMC.h"
36 ClassImp(AliGenReaderCwn);
37
38
39 AliGenReaderCwn::AliGenReaderCwn() 
40 {
41 // Default constructor
42     fNcurrent   = 0;
43     fTreeNtuple = 0;
44 }
45
46 AliGenReaderCwn::~AliGenReaderCwn()
47 {
48     delete fTreeNtuple;
49 }
50
51 void AliGenReaderCwn::Init() 
52 {
53 //
54 // reset the existing file environment and open a new root file if
55 // the pointer to the Fluka tree is null
56     
57     TFile *pFile=0;
58     if (!pFile) {
59         pFile = new TFile(fFileName);
60         pFile->cd();
61         printf("\n I have opened %s file \n", fFileName);
62     }
63 // get the tree address in the Fluka boundary source file
64     fTreeNtuple = (TTree*)gDirectory->Get("h888");
65
66     TTree *h2=fTreeNtuple;
67 //Set branch addresses
68     h2->SetBranchAddress("Nihead",&fNihead);
69     h2->SetBranchAddress("Ihead",fIhead);
70     h2->SetBranchAddress("Nrhead",&fNrhead);
71     h2->SetBranchAddress("Rhead",fRhead);
72     h2->SetBranchAddress("Idpart",&fIdpart);
73     h2->SetBranchAddress("Theta",&fTheta);
74     h2->SetBranchAddress("Phi",&fPhi);
75     h2->SetBranchAddress("P",&fP);
76     h2->SetBranchAddress("E",&fE);
77 }
78
79 Int_t AliGenReaderCwn::NextEvent() 
80 {
81 // Read the next event  
82     Int_t nTracks;
83     fNparticle = 0;
84     TFile* pFile = fTreeNtuple->GetCurrentFile();
85     pFile->cd();
86
87     Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
88     if (fNcurrent < nentries) {
89         Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent);
90         fNcurrent++;
91         
92         Int_t i5=fIhead[4];
93         Int_t i6=fIhead[5];
94         if (i5==0) {
95             printf("\n This should never happen !\n");
96             nTracks = 0;
97         } else {
98             printf("\n Next event contains %d tracks! \n", i6);
99             nTracks = i6;
100         }    
101         fNparticleMax = nTracks;
102         return nTracks;
103     } else {
104         return 0;
105     }
106     return 0;
107 }
108
109 TParticle* AliGenReaderCwn::NextParticle() 
110 {
111 //
112 //  
113     Float_t prwn;
114     Float_t p[4];
115 // Read the next particle
116     if (fCode == kGEANT3) fIdpart=gMC->PDGFromId(fIdpart);
117     Double_t amass = TDatabasePDG::Instance()->GetParticle(fIdpart)->Mass();
118     if(fE<=amass) {
119         Warning("Generate","Particle %d  E = %f mass = %f %f %f \n",
120                 fIdpart,fE,amass, fPhi, fTheta);
121         prwn=0;
122     } else {
123         prwn=sqrt((fE+amass)*(fE-amass));
124     }
125
126     fTheta *= TMath::Pi()/180.;
127     fPhi    = (fPhi-180)*TMath::Pi()/180.;      
128     p[0] = prwn*TMath::Sin(fTheta)*TMath::Cos(fPhi);
129     p[1] = prwn*TMath::Sin(fTheta)*TMath::Sin(fPhi);      
130     p[2] = prwn*TMath::Cos(fTheta);
131     p[3] = fE;
132     TParticle* particle = new TParticle(fIdpart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 
133                                         0., 0., 0., 0.);
134     Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent);
135     fNcurrent++;
136     fNparticle++;
137     return particle;
138 }
139
140
141
142 AliGenReaderCwn& AliGenReaderCwn::operator=(const  AliGenReaderCwn& rhs)
143 {
144 // Assignment operator
145     return *this;
146 }
147
148
149