]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenFunction.cxx
c4f8913127169777f7884fd568c81388a1d4b63a
[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
59 #include "AliRun.h"
60 #include "AliESDtrack.h"
61 #include "AliESDVertex.h"
62 #include "AliGenFunction.h"
63
64 ClassImp(AliGenFunction)
65
66 //-----------------------------------------------------------------------------
67 AliGenFunction::AliGenFunction():
68   AliGenerator(),
69   fBkG(0.),  
70   fFMomentum(0),           // momentum distribution function 
71   fFPhi(0),                // phi distribution function
72   fFTheta(0),              // theta distribution function
73   fFPosition(0),           // position distribution function 
74   fFPdg(0),                // pdg distribution function  
75   fRefRadius(0),           // reference radius to be crossed
76   fZmin(0),                // minimal z at reference radius
77   fZmax(0),                // z at reference radius
78   fMaxTrial(10000)         // maximal number of attempts
79 {
80   //
81   // Default constructor
82   //
83   SetNumberParticles(1);
84 }
85 //-----------------------------------------------------------------------------
86 void AliGenFunction::Generate()
87 {
88   //
89   // Generate one muon
90   //
91   Int_t naccepted =0;
92   
93   for (Int_t ipart=0; ipart<fMaxTrial && naccepted<fNpart; ipart++){
94     //
95     //
96     //
97     Float_t mom[3];
98     Float_t posf[3];
99     Double_t pos[3];
100     Int_t pdg;
101     Double_t ptot,pt,  phi, theta; 
102     //
103     ptot     = fFMomentum->GetRandom();
104     phi      = fFPhi->GetRandom();
105     theta    = fFTheta->GetRandom();
106     pt     = ptot*TMath::Sin(theta);
107     mom[0] = pt*TMath::Cos(phi); 
108     mom[1] = pt*TMath::Sin(phi); 
109     mom[2] = ptot*TMath::Cos(theta);
110
111     //
112     fFPosition->GetRandom3(pos[0],pos[1],pos[2]);
113     posf[0]=pos[0];
114     posf[1]=pos[1];
115     posf[2]=pos[2];
116     pdg = TMath::Nint(fFPdg->GetRandom());
117     if (!IntersectCylinder(fRefRadius,fZmin, fZmax, pdg, posf, mom)) continue;
118     //
119     //
120     Float_t polarization[3]= {0,0,0};
121     Int_t nt;
122     PushTrack(fTrackIt,-1,pdg,mom, posf, polarization,0,kPPrimary,nt);
123     naccepted++;
124   }
125
126
127   return;
128 }
129 //-----------------------------------------------------------------------------
130 void AliGenFunction::Init()
131 {
132   // 
133   // Initialisation, check consistency of selected ranges
134   //
135   printf("************ AliGenFunction ****************\n");
136   printf("************************************************\n");
137   if (!fFMomentum){
138     AliFatal("Momentum distribution function not specified");
139   }
140   if (!fFPosition){
141     AliFatal("Position distribution function not specified");
142   }
143   if (!fFPdg){
144     AliFatal("PDG distribution function not specified");
145   }
146   if (fZmin==fZmax){
147     AliFatal("Z range not specified");
148   }
149
150   return;
151 }
152
153 void AliGenFunction::SetFunctions(TF1 * momentum, TF1 *fphi, TF1 *ftheta,TF3 * position, TF1* pdg){
154   //
155   //
156   //
157   fFMomentum = momentum;
158   fFPhi = fphi;
159   fFTheta = ftheta;
160   fFPosition = position;
161   fFPdg = pdg;
162 }
163
164 void AliGenFunction::SetCylinder(Double_t refR, Double_t zmin, Double_t zmax){
165   //
166   //
167   //
168   fRefRadius = refR;          // reference radius to be crossed
169   fZmin = zmin;               // minimal z at reference radius
170   fZmax = zmax;               // maximal z at reference radius
171
172 }
173
174
175 //-----------------------------------------------------------------------------
176 Bool_t AliGenFunction::IntersectCylinder(Float_t r,Float_t zmin, Float_t zmax,Int_t pdg,
177                                          Float_t o[3],Float_t p[3]) const
178 {
179   //
180   // Intersection between muon and cylinder [-z,+z] with radius r
181   //
182   Double_t mass=0;
183   if (TDatabasePDG::Instance()->GetParticle(pdg)){
184     mass = TDatabasePDG::Instance()->GetParticle(13)->Mass();
185   }
186
187   Float_t en = TMath::Sqrt(mass*mass+p[0]*p[0]+p[1]*p[1]+p[2]*p[2]);
188   TParticle part(pdg,0,0,0,0,0,p[0],p[1],p[2],en,o[0],o[1],o[2],0);
189   AliESDtrack track(&part);
190   Double_t pos[3]={0.,0.,0.},sigma[3]={0.,0.,0.};
191   AliESDVertex origin(pos,sigma);
192
193   track.RelateToVertex(&origin,fBkG,10000.);
194
195   Float_t d0z0[2],covd0z0[3];
196   track.GetImpactParameters(d0z0,covd0z0);
197
198   // check rphi 
199   if(TMath::Abs(d0z0[0])>r) return kFALSE;
200   // check z
201   if(d0z0[1]>zmax) return kFALSE;
202   if(d0z0[1]<zmin) return kFALSE;
203   //
204   return kTRUE;
205 }