]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenExtFile.cxx
Corrected pragma: no need for newObj to access the target data members
[u/mrichter/AliRoot.git] / EVGEN / AliGenExtFile.cxx
CommitLineData
4c039060 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
25454878 16/* $Id$ */
675e9664 17
a48e3376 18// Event generator that using an instance of type AliGenReader
25454878 19// reads particles from a file and applies cuts.
20// Example: In your Config.C you can include the following lines
4f23c932 21// AliGenExtFile *gener = new AliGenExtFile(-1);
25454878 22// gener->SetMomentumRange(0,999);
23// gener->SetPhiRange(-180.,180.);
24// gener->SetThetaRange(0,180);
25// gener->SetYRange(-999,999);
26// AliGenReaderTreeK * reader = new AliGenReaderTreeK();
27// reader->SetFileName("myFileWithTreeK.root");
28// gener->SetReader(reader);
29// gener->Init();
30
675e9664 31
159ec319 32#include <Riostream.h>
084c1b4a 33
3ff7721b 34#include "AliLog.h"
693caace 35#include "AliGenExtFile.h"
4f23c932 36#include "AliRunLoader.h"
9ad9d6f5 37#include "AliHeader.h"
aa4e331d 38#include "AliStack.h"
9ad9d6f5 39#include "AliGenEventHeader.h"
4a33c50d 40#include "AliGenReader.h"
5c3fd7ea 41
a48e3376 42#include <TParticle.h>
693caace 43#include <TFile.h>
a48e3376 44#include <TTree.h>
45
5c3fd7ea 46
83cb58d0 47using std::cout;
48using std::endl;
198bb1c7 49ClassImp(AliGenExtFile)
380956c5 50
51AliGenExtFile::AliGenExtFile()
1c56e311 52 :AliGenMC(),
53 fFileName(0),
9d188d33 54 fReader(0),
55 fStartEvent(0)
693caace 56{
f87cfe57 57// Constructor
693caace 58//
59// Read all particles
e1b35da7 60 fNpart = -1;
693caace 61}
62
63AliGenExtFile::AliGenExtFile(Int_t npart)
1c56e311 64 :AliGenMC(npart),
65 fFileName(0),
9d188d33 66 fReader(0),
67 fStartEvent(0)
693caace 68{
f87cfe57 69// Constructor
a48e3376 70 fName = "ExtFile";
71 fTitle = "Primaries from ext. File";
693caace 72}
73
74//____________________________________________________________
75AliGenExtFile::~AliGenExtFile()
76{
f87cfe57 77// Destructor
a48e3376 78 delete fReader;
693caace 79}
80
a48e3376 81//___________________________________________________________
82void AliGenExtFile::Init()
693caace 83{
a48e3376 84// Initialize
85 if (fReader) fReader->Init();
693caace 86}
87
9d188d33 88//___________________________________________________________
693caace 89void AliGenExtFile::Generate()
90{
f87cfe57 91// Generate particles
693caace 92
eed0bb9a 93 Double_t polar[3] = {0,0,0};
693caace 94 //
eed0bb9a 95 Double_t origin[3] = {0,0,0};
21391258 96 Double_t time = 0.;
1628e931 97 Double_t p[4];
693caace 98 Float_t random[6];
88cb7938 99 Int_t i = 0, j, nt;
693caace 100 //
4f23c932 101 //
102 if (fVertexSmear == kPerEvent) Vertex();
693caace 103
9d188d33 104 // Fast forward up to start Event
105 for (Int_t ie=0; ie<fStartEvent; ++ie ) {
106 Int_t nTracks = fReader->NextEvent();
107 if (nTracks == 0) {
108 // printf("\n No more events !!! !\n");
109 Warning("AliGenExtFile::Generate","\nNo more events in external file!!!\nLast event may be empty or incomplete.\n");
110 return;
111 }
3ff7721b 112 for (Int_t i=0; i<nTracks; ++i) {
113 if (NULL == fReader->NextParticle())
114 AliFatal("Error while skipping tracks");
115 }
9d188d33 116 cout << "Skipping event " << ie << endl;
117 }
3ff7721b 118 fStartEvent = 0; // do not skip events the second time
9d188d33 119
380956c5 120 while(1) {
121 Int_t nTracks = fReader->NextEvent();
122 if (nTracks == 0) {
123 // printf("\n No more events !!! !\n");
124 Warning("AliGenExtFile::Generate","\nNo more events in external file!!!\nLast event may be empty or incomplete.\n");
693caace 125 return;
380956c5 126 }
127
128 //
129 // Particle selection loop
130 //
4f23c932 131 // The selection criterium for the external file generator is as follows:
380956c5 132 //
4f23c932 133 // 1) All tracks are subject to the cuts defined by AliGenerator, i.e.
380956c5 134 // fThetaMin, fThetaMax, fPhiMin, fPhiMax, fPMin, fPMax, fPtMin, fPtMax,
135 // fYMin, fYMax.
136 // If the particle does not satisfy these cuts, it is not put on the
137 // stack.
138 // 2) If fCutOnChild and some specific child is selected (e.g. if
139 // fForceDecay==kSemiElectronic) the event is rejected if NOT EVEN ONE
140 // child falls into the child-cuts.
88cb7938 141 TParticle* iparticle = 0x0;
142
380956c5 143 if(fCutOnChild) {
144 // Count the selected children
145 Int_t nSelected = 0;
4f23c932 146 while ((iparticle=fReader->NextParticle()) ) {
147 Int_t kf = CheckPDGCode(iparticle->GetPdgCode());
148 kf = TMath::Abs(kf);
149 if (ChildSelected(kf) && KinematicSelection(iparticle, 1)) {
150 nSelected++;
151 }
380956c5 152 }
153 if (!nSelected) continue; // No particle selected: Go to next event
154 fReader->RewindEvent();
155 }
a48e3376 156
0d1e7588 157 //
158 // Stack selection loop
159 //
160 class SelectorLogic { // need to do recursive back tracking, requires a "nested" function
161 private:
5263658e 162 Int_t idCount;
0d1e7588 163 map<Int_t,Int_t> firstMotherMap;
164 map<Int_t,Int_t> secondMotherMap;
165 map<Int_t,Bool_t> selectedIdMap;
5263658e 166 map<Int_t,Int_t> newIdMap;
0d1e7588 167 void selectMothersToo(Int_t particleId) {
168 Int_t mum1 = firstMotherMap[particleId];
169 if (mum1 > -1 && !selectedIdMap[mum1]) {
170 selectedIdMap[mum1] = true;
171 selectMothersToo(mum1);
172 }
173 Int_t mum2 = secondMotherMap[particleId];
174 if (mum2 > -1 && !selectedIdMap[mum2]) {
175 selectedIdMap[mum2] = true;
176 selectMothersToo(mum2);
177 }
178 }
179 public:
5263658e 180 void init() {
181 idCount = 0;
182 }
0d1e7588 183 void setData(Int_t id, Int_t mum1, Int_t mum2, Bool_t selected) {
5263658e 184 idCount++; // we know that this function is called in succession of ids, so counting is fine to determine max id
0d1e7588 185 firstMotherMap[id] = mum1;
186 secondMotherMap[id] = mum2;
187 selectedIdMap[id] = selected;
5263658e 188 }
189 void reselectCuttedMothersAndRemapIDs() {
190 for (Int_t id = 0; id < idCount; ++id) {
191 if (selectedIdMap[id]) {
192 selectMothersToo(id);
193 }
194 }
195 Int_t newId = 0;
196 for (Int_t id = 0; id < idCount; id++) {
197 if (selectedIdMap[id]) {
198 newIdMap[id] = newId; ++newId;
396b3bc9 199 } else {
200 newIdMap[id] = -1;
5263658e 201 }
0d1e7588 202 }
203 }
204 Bool_t isSelected(Int_t id) {
205 return selectedIdMap[id];
206 }
5263658e 207 Int_t newId(Int_t id) {
208 if (id == -1) return -1;
209 return newIdMap[id];
210 }
0d1e7588 211 };
212 SelectorLogic selector;
5263658e 213 selector.init();
0d1e7588 214 for (Int_t i = 0; i < nTracks; i++) {
215 TParticle* jparticle = fReader->NextParticle();
216 selector.setData(i,
217 jparticle->GetFirstMother(),
218 jparticle->GetSecondMother(),
219 KinematicSelection(jparticle,0));
220 }
5263658e 221 selector.reselectCuttedMothersAndRemapIDs();
0d1e7588 222 fReader->RewindEvent();
223
380956c5 224 //
225 // Stack filling loop
226 //
7cf36cae 227 fNprimaries = 0;
380956c5 228 for (i = 0; i < nTracks; i++) {
bce6240d 229 TParticle* jparticle = fReader->NextParticle();
0d1e7588 230 Bool_t selected = selector.isSelected(i);
bce6240d 231 if (!selected) {
bce6240d 232 continue;
233 }
5263658e 234 Int_t parent = selector.newId(jparticle->GetFirstMother());
235// printf("particle %d -> %d, with mother %d -> %d\n", i, selector.newId(i), jparticle->GetFirstMother(), parent);
bce6240d 236
ea79897e 237 p[0] = jparticle->Px();
238 p[1] = jparticle->Py();
239 p[2] = jparticle->Pz();
1628e931 240 p[3] = jparticle->Energy();
241
b7a824ae 242 Int_t idpart = jparticle->GetPdgCode();
4f23c932 243 if(fVertexSmear==kPerTrack)
244 {
245 Rndm(random,6);
246 for (j = 0; j < 3; j++) {
247 origin[j]=fOrigin[j]+
248 fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
249 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
250 }
21391258 251 Rndm(random,2);
252 time = fTimeOrigin + fOsigma[2]/TMath::Ccgs()*
253 TMath::Cos(2*random[0]*TMath::Pi())*
254 TMath::Sqrt(-2*TMath::Log(random[1]));
380956c5 255 } else {
ea79897e 256 origin[0] = fVertex[0] + jparticle->Vx();
257 origin[1] = fVertex[1] + jparticle->Vy();
258 origin[2] = fVertex[2] + jparticle->Vz();
21391258 259 time = fTime + jparticle->T();
380956c5 260 }
8f126362 261 Int_t doTracking = fTrackIt && selected && (jparticle->TestBit(kTransportBit));
4f23c932 262
eed0bb9a 263 PushTrack(doTracking, parent, idpart,
21391258 264 p[0], p[1], p[2], p[3], origin[0], origin[1], origin[2], time,
eed0bb9a 265 polar[0], polar[1], polar[2],
266 kPPrimary, nt, 1., jparticle->GetStatusCode());
267
4f23c932 268 KeepTrack(nt);
7cf36cae 269 fNprimaries++;
380956c5 270 } // track loop
0b05fbf9 271
272 // Generated event header
7a005d69 273 AliGenEventHeader * header = fReader->GetGenEventHeader();
274 if (!header) header = new AliGenEventHeader();
7cf36cae 275 header->SetNProduced(fNprimaries);
0b05fbf9 276 header->SetPrimaryVertex(fVertex);
21391258 277 header->SetInteractionTime(fTime);
7cf36cae 278 AddHeader(header);
380956c5 279 break;
4f23c932 280
380956c5 281 } // event loop
4f23c932 282
380956c5 283 SetHighWaterMark(nt);
2c354237 284 CdEventFile();
285}
380956c5 286
9d188d33 287//___________________________________________________________
2c354237 288void AliGenExtFile::CdEventFile()
289{
290// CD back to the event file
33c3c91a 291 AliRunLoader::Instance()->CdGAFile();
693caace 292}
293
294
693caace 295
296
297