]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenExtFile.cxx
Typo corrected.
[u/mrichter/AliRoot.git] / EVGEN / AliGenExtFile.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 $Log$
18 Revision 1.21  2002/04/26 10:39:31  morsch
19 AliGenExtFile derives from AliGenMC. Generate() uses methods from AliGenMC (N. Carrer)
20
21 Revision 1.20  2002/03/22 09:43:28  morsch
22 Don't delete the TParticle.
23
24 Revision 1.19  2002/02/08 16:50:50  morsch
25 Add name and title in constructor.
26
27 Revision 1.18  2001/11/12 14:31:00  morsch
28 Memory leaks fixed. (M. Bondila)
29
30 Revision 1.17  2001/11/09 09:12:58  morsch
31 Generalization by using AliGenReader object to read particles from file.
32
33 Revision 1.16  2001/07/27 17:09:35  morsch
34 Use local SetTrack, KeepTrack and SetHighWaterMark methods
35 to delegate either to local stack or to stack owned by AliRun.
36 (Piotr Skowronski, A.M.)
37
38 Revision 1.15  2001/01/23 13:29:37  morsch
39 Add method SetParticleCode and enum type Code_t to handle both PDG (new ntuples)
40 and GEANT3 codes (old ntuples) in input file.
41
42 Revision 1.14  2000/12/21 16:24:06  morsch
43 Coding convention clean-up
44
45 Revision 1.13  2000/11/30 07:12:50  alibrary
46 Introducing new Rndm and QA classes
47
48 Revision 1.12  2000/10/27 13:54:45  morsch
49 Remove explicite reference to input file from constuctor.
50
51 Revision 1.11  2000/10/02 21:28:06  fca
52 Removal of useless dependecies via forward declarations
53
54 Revision 1.10  2000/07/11 18:24:55  fca
55 Coding convention corrections + few minor bug fixes
56
57 Revision 1.9  2000/06/14 15:20:09  morsch
58 Include clean-up (IH)
59
60 Revision 1.8  2000/06/09 20:36:44  morsch
61 All coding rule violations except RS3 corrected
62
63 Revision 1.7  2000/02/16 14:56:27  morsch
64 Convert geant particle code into pdg code before putting particle on the stack.
65
66 Revision 1.6  1999/11/09 07:38:48  fca
67 Changes for compatibility with version 2.23 of ROOT
68
69 Revision 1.5  1999/09/29 09:24:12  fca
70 Introduction of the Copyright and cvs Log
71
72 */
73
74
75 // Event generator that using an instance of type AliGenReader
76 // reads particles from a file and applies cuts. 
77
78 #include <iostream.h>
79
80 #include "AliGenExtFile.h"
81 #include "AliRun.h"
82
83 #include <TParticle.h>
84 #include <TFile.h>
85 #include <TTree.h>
86
87
88  ClassImp(AliGenExtFile)
89
90 AliGenExtFile::AliGenExtFile()
91   :AliGenMC()
92 {
93 //  Constructor
94 //
95 //  Read all particles
96     fNpart  =- 1;
97     fReader =  0;
98 }
99
100 AliGenExtFile::AliGenExtFile(Int_t npart)
101     :AliGenMC(npart)
102 {
103 //  Constructor
104     fName   = "ExtFile";
105     fTitle  = "Primaries from ext. File";
106     fReader = 0;
107 }
108
109 AliGenExtFile::AliGenExtFile(const AliGenExtFile & ExtFile)
110 {
111 // copy constructor
112 }
113 //____________________________________________________________
114 AliGenExtFile::~AliGenExtFile()
115 {
116 // Destructor
117     delete fReader;
118 }
119
120 //___________________________________________________________
121 void AliGenExtFile::Init()
122 {
123 // Initialize
124     if (fReader) fReader->Init();
125 }
126
127     
128 void AliGenExtFile::Generate()
129 {
130 // Generate particles
131
132   Float_t polar[3]  = {0,0,0};
133   //
134   Float_t origin[3] = {0,0,0};
135   Float_t p[3];
136   Float_t random[6];
137   Int_t i, j, nt;
138   //
139   for (j=0;j<3;j++) origin[j]=fOrigin[j];
140   if(fVertexSmear == kPerTrack) {
141     Rndm(random,6);
142     for (j = 0; j < 3; j++) {
143         origin[j] += fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
144             TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
145     }
146   }
147
148   while(1) {
149     Int_t nTracks = fReader->NextEvent();       
150     if (nTracks == 0) {
151       // printf("\n No more events !!! !\n");
152       Warning("AliGenExtFile::Generate","\nNo more events in external file!!!\nLast event may be empty or incomplete.\n");
153       return;
154     }
155
156     //
157     // Particle selection loop
158     //
159     // The selction criterium for the external file generator is as follows:
160     //
161     // 1) All tracs are subjects to the cuts defined by AliGenerator, i.e.
162     //    fThetaMin, fThetaMax, fPhiMin, fPhiMax, fPMin, fPMax, fPtMin, fPtMax,
163     //    fYMin, fYMax.
164     //    If the particle does not satisfy these cuts, it is not put on the
165     //    stack.
166     // 2) If fCutOnChild and some specific child is selected (e.g. if
167     //    fForceDecay==kSemiElectronic) the event is rejected if NOT EVEN ONE
168     //    child falls into the child-cuts.
169     if(fCutOnChild) {
170       // Count the selected children
171       Int_t nSelected = 0;
172
173       for (i = 0; i < nTracks; i++) {
174         TParticle* iparticle = fReader->NextParticle();
175         Int_t kf = CheckPDGCode(iparticle->GetPdgCode());
176         kf = TMath::Abs(kf);
177         if (ChildSelected(kf) && KinematicSelection(iparticle, 1)) {
178           nSelected++;
179         }
180       }
181       if (!nSelected) continue;    // No particle selected:  Go to next event
182       fReader->RewindEvent();
183     }
184
185     //
186     // Stack filling loop
187     //
188     for (i = 0; i < nTracks; i++) {
189
190       TParticle* iparticle = fReader->NextParticle();
191       if (!KinematicSelection(iparticle,0)) {
192         Double_t  pz   = iparticle->Pz();
193         Double_t  e    = iparticle->Energy();
194         Double_t  y;
195         if ((e-pz) == 0) {
196           y = 20.;
197         } else if ((e+pz) == 0.) {
198           y = -20.;
199         } else {
200           y = 0.5*TMath::Log((e+pz)/(e-pz));      
201         }
202         printf("\n Not selected %d %f %f %f %f %f", i,
203                iparticle->Theta(),
204                iparticle->Phi(),
205                iparticle->P(),
206                iparticle->Pt(),
207                y);
208         delete iparticle;
209         continue;
210       }
211       p[0] = iparticle->Px();
212       p[1] = iparticle->Py();
213       p[2] = iparticle->Pz();
214       Int_t idpart = iparticle->GetPdgCode();
215       if(fVertexSmear==kPerTrack) {
216           Rndm(random,6);
217           for (j = 0; j < 3; j++) {
218               origin[j]=fOrigin[j]
219                   +fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
220                   TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
221           }
222       }
223       Int_t decayed = iparticle->GetFirstDaughter();
224       Int_t doTracking = fTrackIt && (decayed < 0) &&
225                          (TMath::Abs(idpart) > 10);
226       // printf("*** pdg, first daughter, trk = %d, %d, %d\n",
227       //   idpart,decayed, doTracking);
228       SetTrack(doTracking,-1,idpart,p,origin,polar,0,kPPrimary,nt);
229       KeepTrack(nt);
230     } // track loop
231
232     break;
233
234   } // event loop
235
236   SetHighWaterMark(nt);
237   CdEventFile();
238 }
239
240 void AliGenExtFile::CdEventFile()
241 {
242 // CD back to the event file
243   TFile *pFile=0;
244   if (!gAlice) {
245       gAlice = (AliRun*)pFile->Get("gAlice");
246       if (gAlice) printf("AliRun object found on file\n");
247       if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
248   }
249   TTree *fAli=gAlice->TreeK();
250   if (fAli) pFile =fAli->GetCurrentFile();
251   pFile->cd();    
252 }
253
254
255 AliGenExtFile& AliGenExtFile::operator=(const  AliGenExtFile& rhs)
256 {
257 // Assignment operator
258     return *this;
259 }
260
261
262
263
264
265
266
267