]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenSlowNucleons.cxx
Option kBJpsi added.
[u/mrichter/AliRoot.git] / EVGEN / AliGenSlowNucleons.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 /* $Id$ */
17
18 //
19 //  Generator for slow nucleons in pA interactions. 
20 //  Source is modelled by a relativistic Maxwell distributions.
21 //  This class cooparates with AliCollisionGeometry if used inside AliGenCocktail.
22 //  In this case the number of slow nucleons is determined from the number of wounded nuclei
23 //  using a realisation of AliSlowNucleonModel.
24 //  Original code by  Ferenc Sikler  <sikler@rmki.kfki.hu>
25 // 
26
27 #include <TDatabasePDG.h>
28 #include <TPDGCode.h>
29 #include <TH2F.h>
30 #include <TCanvas.h>
31
32 #include "AliCollisionGeometry.h"
33 #include "AliGenSlowNucleons.h"
34 #include "AliSlowNucleonModel.h"
35
36  ClassImp(AliGenSlowNucleons)
37     
38  AliGenSlowNucleons::AliGenSlowNucleons():AliGenerator(-1)
39 {
40 // Default constructor
41     fSlowNucleonModel = 0;
42     fCollisionGeometry = 0;
43 }
44
45 AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
46     :AliGenerator(npart)
47 {
48 // Constructor
49     fName  = "SlowNucleons";
50     fTitle = "Generator for gray particles in pA collisions";
51     SetPmax();
52     SetTarget();
53     SetNominalCmsEnergy();
54     SetCharge();
55     SetTemperature();
56     SetBetaSource();
57     fSlowNucleonModel = new AliSlowNucleonModel();
58     fCollisionGeometry = 0;
59     fDebug = 0;
60 }
61
62 AliGenSlowNucleons::AliGenSlowNucleons(const AliGenSlowNucleons & sn):
63     AliGenerator(sn)
64 {
65 // Copy constructor
66     sn.Copy(*this);
67 }
68
69 //____________________________________________________________
70 AliGenSlowNucleons::~AliGenSlowNucleons()
71 {
72 // Destructor
73     delete  fSlowNucleonModel;
74 }
75
76
77 void AliGenSlowNucleons::Init()
78 {
79   //
80   // Initialization
81   //
82     Float_t kMass  = TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
83     fMomentum = fCMS/2. * fZTarget / fATarget;
84     fBeta     = fMomentum / TMath::Sqrt(kMass * kMass + fMomentum * fMomentum);
85     if (fDebug) {
86         fDebugHist1 = new TH2F("DebugHist1", "nu vs N_slow", 100, 0., 100., 20, 0., 20.);
87         fDebugHist2 = new TH2F("DebugHist2", "b  vs N_slow", 100, 0., 100., 15, 0., 15.);
88     }
89 }
90
91 void AliGenSlowNucleons::FinishRun()
92 {
93 // End of run action
94 // Show histogram for debugging if requested.
95     if (fDebug) {
96         TCanvas *c = new TCanvas("c","Canvas 1",400,10,600,700);
97         c->Divide(2,1);
98         c->cd(1);
99         fDebugHist1->Draw();
100         c->cd(2);
101         fDebugHist2->Draw();
102     }
103 }
104
105
106 void AliGenSlowNucleons::Generate()
107 {
108   //
109   // Generate one event
110   //
111   //
112   // Communication with Gray Particle Model 
113   // 
114     if (fCollisionGeometry) {
115         Float_t b   = fCollisionGeometry->ImpactParameter();
116         Int_t  nn   = fCollisionGeometry->NN();
117         Int_t  nwn  = fCollisionGeometry->NwN();
118         Int_t  nnw  = fCollisionGeometry->NNw();
119         Int_t  nwnw = fCollisionGeometry->NwNw();
120
121         fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn);
122         if (fDebug) {
123             printf("Nucleons %d %d %d %d \n", fNgp, fNgn, fNbp, fNbn);
124             fDebugHist1->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), fCollisionGeometry->NwN(), 1.);
125             fDebugHist2->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), b, 1.);
126             printf("AliGenSlowNucleons: Impact parameter from Collision Geometry %f %d %d %d %d\n", 
127                    b, nn, nwn, nnw, nwnw);
128         }
129     }     
130
131    //
132     Float_t p[3];
133     Float_t origin[3] = {0., 0., 0.};
134     Float_t polar [3] = {0., 0., 0.};    
135     Int_t nt, i, j;
136     Int_t kf;
137     
138     if(fVertexSmear == kPerEvent) {
139         Vertex();
140         for (j=0; j < 3; j++) origin[j] = fVertex[j];
141     } // if kPerEvent
142 //
143 //  Gray protons
144 //
145     fCharge = 1;
146     kf = kProton;    
147     for(i = 0; i < fNgp; i++) {
148         GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
149         PushTrack(fTrackIt, -1, kf, p, origin, polar,
150                  0., kPNoProcess, nt, 1.);
151         KeepTrack(nt);
152     }
153 //
154 //  Gray neutrons
155 //
156     fCharge = 0;
157     kf = kNeutron;    
158     for(i = 0; i < fNgn; i++) {
159         GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
160         PushTrack(fTrackIt, -1, kf, p, origin, polar,
161                  0., kPNoProcess, nt, 1.);
162         KeepTrack(nt);
163     }
164 //
165 //  Black protons
166 //
167     fCharge = 1;
168     kf = kProton;    
169     for(i = 0; i < fNbp; i++) {
170         GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
171         PushTrack(fTrackIt, -1, kf, p, origin, polar,
172                  0., kPNoProcess, nt, 1.);
173         KeepTrack(nt);
174     }
175 //
176 //  Black neutrons
177 //
178     fCharge = 0;
179     kf = kNeutron;    
180     for(i = 0; i < fNbn; i++) {
181         GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
182         PushTrack(fTrackIt, -1, kf, p, origin, polar,
183                  0., kPNoProcess, nt, 1.);
184         KeepTrack(nt);
185     }
186 }
187
188
189
190
191 void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, Double_t beta, Float_t* q)
192
193 {
194 /* 
195    Emit a slow nucleon with "temperature" T [GeV], 
196    from a source moving with velocity beta         
197    Three-momentum [GeV/c] is given back in q[3]    
198 */
199
200  Double_t m, pmax, p, f, theta, phi;
201  TDatabasePDG * pdg = TDatabasePDG::Instance();
202  const Double_t kMassProton  = pdg->GetParticle(kProton) ->Mass();
203  const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
204  
205  /* Select nucleon type */
206  if(charge == 0) m = kMassNeutron;
207  else m = kMassProton;
208
209  /* Momentum at maximum of Maxwell-distribution */
210
211  pmax = TMath::Sqrt(2*T*(T+sqrt(T*T+m*m)));
212
213  /* Try until proper momentum                                  */
214  /* for lack of primitive function of the Maxwell-distribution */
215  /* a brute force trial-accept loop, normalized at pmax        */
216
217  do
218  {
219      p = Rndm() * fPmax;
220      f = Maxwell(m, p, T) / Maxwell(m , pmax, T);
221  }
222  while(f < Rndm());
223
224  /* Spherical symmetric emission */
225  theta = TMath::ACos(2. * Rndm() - 1.);
226  phi   = 2. * TMath::Pi() * Rndm();
227
228  /* Determine momentum components in system of the moving source */
229  q[0] = p * TMath::Sin(theta) * TMath::Cos(phi);
230  q[1] = p * TMath::Sin(theta) * TMath::Sin(phi);
231  q[2] = p * TMath::Cos(theta);
232
233  /* Transform to system of the target nucleus                             */
234  /* beta is passed as negative, because the gray nucleons are slowed down */
235  Lorentz(m, -beta, q);
236
237  /* Transform to laboratory system */
238  Lorentz(m, fBeta, q);
239 }
240
241 Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T)
242 {
243 /* Relativistic Maxwell-distribution */
244     Double_t ekin;
245     ekin = TMath::Sqrt(p*p+m*m)-m;
246     return (p*p * exp(-ekin/T));
247 }
248
249
250 void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q)
251 {
252 /* Lorentz transform in the direction of q[2] */
253  
254     Double_t gamma  = 1/sqrt(1-beta*beta); 
255     Double_t energy = sqrt(m*m + q[0]*q[0] + q[1]*q[1] + q[2]*q[2]);
256     q[2] = gamma * (q[2] + beta*energy);
257 }
258
259           
260 AliGenSlowNucleons& AliGenSlowNucleons::operator=(const  AliGenSlowNucleons& rhs)
261 {
262 // Assignment operator
263     rhs.Copy(*this);
264     return *this;
265 }
266
267 void AliGenSlowNucleons::Copy(TObject&) const
268 {
269     //
270     // Copy 
271     //
272     Fatal("Copy","Not implemented!\n");
273 }
274
275
276
277
278
279
280