]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenFunction.cxx
Fix for bug #75588
[u/mrichter/AliRoot.git] / EVGEN / AliGenFunction.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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 // Generator for particles according generic functions
17 //  
18 //  TF1 *   fFMomentum;           // momentum distribution function inGeV
19 //  TF1 *   fFPhi;                // phi distribution function in rad
20 //  TF1 *   fFTheta;              // theta distribution function in rad
21 //  TF3 *   fFPosition;           // position distribution function in cm
22 //  TF1 *   fFPdg;                // pdg distribution function  
23 //  We assume that the moment, postion and PDG code of particles are independent  
24 //  Only tracks/particle crossing the reference radius at given z range
25 //
26 // Origin: marian.ivanov@cern.ch
27
28
29 /*
30   Example generic function for cosmic generation:
31  
32 //
33         AliGenFunction *generCosmic = new AliGenFunction;
34         generCosmic->SetNumberParticles(100);
35         TF1 *fmom   = new TF1("fmom","TMath::Exp(-x/8)",0,30);
36         //
37         TF1 *fphi   = new TF1("fphi","TMath::Gaus(x,-TMath::Pi()/2,0.3)",-3.14,3.14);
38         TF1 *ftheta = new TF1("ftheta","TMath::Gaus(x,TMath::Pi()/2,0.3)",-3.14,3.14);
39         TF3 *fpos = new TF3("fpos","1+(x+y+z)*0",-2000,2000,700,701,-2000,2000);
40         TF1 * fpdg= new TF1("fpdg","1*(abs(x-13)<0.1)+1*(abs(x+13)<0.1)",-300,300);
41         fpdg->SetNpx(10000);  // neccessary - number of bins for generation
42         generCosmic->SetVertexSmear(kPerEvent);
43         generCosmic->SetFunctions(fmom,fphi,ftheta,fpos,fpdg);
44         generCosmic->SetCylinder(375,-375,375);
45         generCosmic->SetBkG(0.2);
46         gener->AddGenerator(generCosmic,"generCosmic",1);
47
48
49   
50 */
51
52
53
54
55 #include <TParticle.h>
56 #include <TF1.h>
57 #include <TF3.h>
58 #include <TDatabasePDG.h>
59
60 #include "AliRun.h"
61 #include "AliESDtrack.h"
62 #include "AliESDVertex.h"
63 #include "AliGenFunction.h"
64 #include "AliGenEventHeader.h"
65
66 ClassImp(AliGenFunction)
67
68 //-----------------------------------------------------------------------------
69 AliGenFunction::AliGenFunction():
70   AliGenerator(),
71   fBkG(0.),  
72   fFMomentum(0),           // momentum distribution function 
73   fFPhi(0),                // phi distribution function
74   fFTheta(0),              // theta distribution function
75   fFPosition(0),           // position distribution function 
76   fFPdg(0),                // pdg distribution function  
77   fRefRadius(0),           // reference radius to be crossed
78   fZmin(0),                // minimal z at reference radius
79   fZmax(0),                // z at reference radius
80   fMaxTrial(10000)         // maximal number of attempts
81 {
82   //
83   // Default constructor
84   //
85   SetNumberParticles(1);
86 }
87
88 AliGenFunction::AliGenFunction(const AliGenFunction& func):
89   AliGenerator(),
90   fBkG(func.fBkG),  
91   fFMomentum(func.fFMomentum),           // momentum distribution function 
92   fFPhi(func.fFPhi),                     // phi distribution function
93   fFTheta(func.fFTheta),                 // theta distribution function
94   fFPosition(func.fFPosition),           // position distribution function 
95   fFPdg(func.fFPdg),                     // pdg distribution function  
96   fRefRadius(func.fRefRadius),           // reference radius to be crossed
97   fZmin(func.fZmin),                     // minimal z at reference radius
98   fZmax(func.fZmax),                     // z at reference radius
99   fMaxTrial(10000)                       // maximal number of attempts
100 {
101     // Copy constructor
102     SetNumberParticles(1);
103 }
104
105 AliGenFunction & AliGenFunction::operator=(const AliGenFunction& func)
106 {
107     // Assigment operator
108       if(&func == this) return *this;
109       fBkG       = func.fBkG;
110       fFMomentum = func.fFMomentum; 
111       fFPhi      = func.fFPhi;      
112       fFTheta    = func.fFTheta;    
113       fFPosition = func.fFPosition; 
114       fFPdg      = func.fFPdg;      
115       fRefRadius = func.fRefRadius;      
116       fZmin      = func.fZmin;                
117       fZmax      = func.fZmax;                
118       fMaxTrial  = func.fMaxTrial;
119       return *this;
120 }
121
122
123 //-----------------------------------------------------------------------------
124 void AliGenFunction::Generate()
125 {
126   //
127   // Generate one muon
128   //
129   Int_t naccepted =0;
130   
131   for (Int_t ipart=0; ipart<fMaxTrial && naccepted<fNpart; ipart++){
132     //
133     //
134     //
135     Float_t mom[3];
136     Float_t posf[3];
137     Double_t pos[3];
138     Int_t pdg;
139     Double_t ptot, pt,  phi, theta; 
140     //
141     ptot     = fFMomentum->GetRandom();
142     phi      = fFPhi->GetRandom();
143     theta    = fFTheta->GetRandom();
144     pt     = ptot*TMath::Sin(theta);
145     mom[0] = pt*TMath::Cos(phi); 
146     mom[1] = pt*TMath::Sin(phi); 
147     mom[2] = ptot*TMath::Cos(theta);
148
149     //
150     fFPosition->GetRandom3(pos[0],pos[1],pos[2]);
151     posf[0]=pos[0];
152     posf[1]=pos[1];
153     posf[2]=pos[2];
154     pdg = TMath::Nint(fFPdg->GetRandom());
155     if (!IntersectCylinder(fRefRadius,fZmin, fZmax, pdg, posf, mom)) continue;
156     //
157     //
158     Float_t polarization[3]= {0,0,0};
159     Int_t nt;
160     PushTrack(fTrackIt,-1,pdg,mom, posf, polarization,0,kPPrimary,nt);
161     naccepted++;
162   }
163
164   AliGenEventHeader* header = new AliGenEventHeader("THn");
165   gAlice->SetGenEventHeader(header);
166   return;
167 }
168 //-----------------------------------------------------------------------------
169 void AliGenFunction::Init()
170 {
171   // 
172   // Initialisation, check consistency of selected ranges
173   //
174   printf("************ AliGenFunction ****************\n");
175   printf("************************************************\n");
176   if (!fFMomentum){
177     AliFatal("Momentum distribution function not specified");
178   }
179   if (!fFPosition){
180     AliFatal("Position distribution function not specified");
181   }
182   if (!fFPdg){
183     AliFatal("PDG distribution function not specified");
184   }
185   if (fZmin==fZmax){
186     AliFatal("Z range not specified");
187   }
188
189   return;
190 }
191
192 void AliGenFunction::SetFunctions(TF1 * momentum, TF1 *fphi, TF1 *ftheta,TF3 * position, TF1* pdg){
193   //
194   //
195   //
196   fFMomentum = momentum;
197   fFPhi = fphi;
198   fFTheta = ftheta;
199   fFPosition = position;
200   fFPdg = pdg;
201 }
202
203 void AliGenFunction::SetCylinder(Double_t refR, Double_t zmin, Double_t zmax){
204   //
205   //
206   //
207   fRefRadius = refR;          // reference radius to be crossed
208   fZmin = zmin;               // minimal z at reference radius
209   fZmax = zmax;               // maximal z at reference radius
210
211 }
212
213
214 //-----------------------------------------------------------------------------
215 Bool_t AliGenFunction::IntersectCylinder(Float_t r,Float_t zmin, Float_t zmax,Int_t pdg,
216                                          Float_t o[3],Float_t p[3]) const
217 {
218   //
219   // Intersection between muon and cylinder [-z,+z] with radius r
220   //
221   Double_t mass=0;
222   if (TDatabasePDG::Instance()->GetParticle(pdg)){
223     mass = TDatabasePDG::Instance()->GetParticle(13)->Mass();
224   }
225
226   Float_t en = TMath::Sqrt(mass*mass+p[0]*p[0]+p[1]*p[1]+p[2]*p[2]);
227   TParticle part(pdg,0,0,0,0,0,p[0],p[1],p[2],en,o[0],o[1],o[2],0);
228   AliESDtrack track(&part);
229   Double_t pos[3]={0.,0.,0.},sigma[3]={0.,0.,0.};
230   AliESDVertex origin(pos,sigma);
231
232   track.RelateToVertex(&origin,fBkG,10000.);
233
234   Float_t d0z0[2],covd0z0[3];
235   track.GetImpactParameters(d0z0,covd0z0);
236
237   // check rphi 
238   if(TMath::Abs(d0z0[0])>r) return kFALSE;
239   // check z
240   if(d0z0[1]>zmax) return kFALSE;
241   if(d0z0[1]<zmin) return kFALSE;
242   //
243   return kTRUE;
244 }