Line data Source code
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 : // Parameterisation of pi and K, eta and pt distributions
19 : // used for the ALICE TDRs.
20 : // eta: according to HIJING (shadowing + quenching)
21 : // pT : according to CDF measurement at 1.8 TeV
22 : // Author: andreas.morsch@cern.ch
23 :
24 :
25 : //Begin_Html
26 : /*
27 : <img src="picts/AliGeneratorClass.gif">
28 : </pre>
29 : <br clear=left>
30 : <font size=+2 color=red>
31 : <p>The responsible person for this module is
32 : <a href="mailto:andreas.morsch@cern.ch">Andreas Morsch</a>.
33 : </font>
34 : <pre>
35 : */
36 : //End_Html
37 : // //
38 : ///////////////////////////////////////////////////////////////////
39 :
40 : #include <TArrayF.h>
41 : #include <TCanvas.h>
42 : #include <TClonesArray.h>
43 : #include <TDatabasePDG.h>
44 : #include <TF1.h>
45 : #include <TH1.h>
46 : #include <TPDGCode.h>
47 : #include <TParticle.h>
48 : #include <TROOT.h>
49 : #include <TVirtualMC.h>
50 :
51 : #include "AliConst.h"
52 : #include "AliDecayer.h"
53 : #include "AliGenEventHeader.h"
54 : #include "AliGenHIJINGpara.h"
55 : #include "AliLog.h"
56 : #include "AliRun.h"
57 :
58 6 : ClassImp(AliGenHIJINGpara)
59 :
60 :
61 :
62 : //_____________________________________________________________________________
63 : static Double_t ptpi(const Double_t *px, const Double_t *)
64 : {
65 : //
66 : // PT-PARAMETERIZATION CDF, PRL 61(88) 1819
67 : // POWER LAW FOR PT > 500 MEV
68 : // MT SCALING BELOW (T=160 MEV)
69 : //
70 : const Double_t kp0 = 1.3;
71 : const Double_t kxn = 8.28;
72 : const Double_t kxlim = 0.5;
73 : const Double_t kt = 0.160;
74 : const Double_t kxmpi = 0.139;
75 : const Double_t kb = 1.;
76 : Double_t y, y1, xmpi2, ynorm, a;
77 0 : Double_t x = *px;
78 : //
79 0 : y1 = TMath::Power(kp0 / (kp0 + kxlim), kxn);
80 : xmpi2 = kxmpi * kxmpi;
81 0 : ynorm = kb * (TMath::Exp(-sqrt(kxlim * kxlim + xmpi2) / kt ));
82 0 : a = ynorm / y1;
83 0 : if (x > kxlim)
84 0 : y = a * TMath::Power(kp0 / (kp0 + x), kxn);
85 : else
86 0 : y = kb* TMath::Exp(-sqrt(x * x + xmpi2) / kt);
87 :
88 0 : return y*x;
89 : }
90 :
91 : //_____________________________________________________________________________
92 : static Double_t ptscal(Double_t pt, Int_t np)
93 : {
94 : // SCALING EN MASSE PAR RAPPORT A PTPI
95 : // MASS PI,K,ETA,RHO,OMEGA,ETA',PHI
96 : const Double_t khm[10] = {.13957,.493,.5488,.769,.7826,.958,1.02,0,0,0};
97 : // VALUE MESON/PI AT 5 GEV
98 : const Double_t kfmax[10]={1.,0.3,0.55,1.0,1.0,1.0,1.0,0,0,0};
99 0 : np--;
100 0 : Double_t f5=TMath::Power(((
101 0 : sqrt(100.018215)+2.)/(sqrt(100.+khm[np]*khm[np])+2.0)),12.3);
102 0 : Double_t fmax2=f5/kfmax[np];
103 : // PIONS
104 0 : Double_t ptpion=100.*ptpi(&pt, (Double_t*) 0);
105 0 : Double_t fmtscal=TMath::Power(((
106 0 : sqrt(pt*pt+0.018215)+2.)/ (sqrt(pt*pt+khm[np]*khm[np])+2.0)),12.3)/
107 : fmax2;
108 0 : return fmtscal*ptpion;
109 : }
110 :
111 : //_____________________________________________________________________________
112 : static Double_t ptka( Double_t *px, Double_t *)
113 : {
114 : //
115 : // pt parametrisation for k
116 : //
117 0 : return ptscal(*px,2);
118 : }
119 :
120 :
121 : //_____________________________________________________________________________
122 : static Double_t etapic( Double_t *py, Double_t *)
123 : {
124 : //
125 : // eta parametrisation for pi
126 : //
127 : const Double_t ka1 = 4913.;
128 : const Double_t ka2 = 1819.;
129 : const Double_t keta1 = 0.22;
130 : const Double_t keta2 = 3.66;
131 : const Double_t kdeta1 = 1.47;
132 : const Double_t kdeta2 = 1.51;
133 0 : Double_t y=TMath::Abs(*py);
134 : //
135 0 : Double_t ex1 = (y-keta1)*(y-keta1)/(2*kdeta1*kdeta1);
136 0 : Double_t ex2 = (y-keta2)*(y-keta2)/(2*kdeta2*kdeta2);
137 0 : return ka1*TMath::Exp(-ex1)+ka2*TMath::Exp(-ex2);
138 : }
139 :
140 : //_____________________________________________________________________________
141 : static Double_t etakac( Double_t *py, Double_t *)
142 : {
143 : //
144 : // eta parametrisation for ka
145 : //
146 : const Double_t ka1 = 497.6;
147 : const Double_t ka2 = 215.6;
148 : const Double_t keta1 = 0.79;
149 : const Double_t keta2 = 4.09;
150 : const Double_t kdeta1 = 1.54;
151 : const Double_t kdeta2 = 1.40;
152 0 : Double_t y=TMath::Abs(*py);
153 : //
154 0 : Double_t ex1 = (y-keta1)*(y-keta1)/(2*kdeta1*kdeta1);
155 0 : Double_t ex2 = (y-keta2)*(y-keta2)/(2*kdeta2*kdeta2);
156 0 : return ka1*TMath::Exp(-ex1)+ka2*TMath::Exp(-ex2);
157 : }
158 :
159 : //_____________________________________________________________________________
160 : AliGenHIJINGpara::AliGenHIJINGpara()
161 0 : :AliGenerator(),
162 0 : fNt(-1),
163 0 : fNpartProd(0),
164 0 : fPi0Decays(kFALSE),
165 0 : fPtWgtPi(0.),
166 0 : fPtWgtKa(0.),
167 0 : fPtpi(0),
168 0 : fPtka(0),
169 0 : fETApic(0),
170 0 : fETAkac(0),
171 0 : fDecayer(0)
172 0 : {
173 : //
174 : // Default constructor
175 : //
176 0 : SetCutVertexZ();
177 0 : SetPtRange();
178 0 : }
179 :
180 : //_____________________________________________________________________________
181 : AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
182 0 : :AliGenerator(npart),
183 0 : fNt(-1),
184 0 : fNpartProd(npart),
185 0 : fPi0Decays(kFALSE),
186 0 : fPtWgtPi(0.),
187 0 : fPtWgtKa(0.),
188 0 : fPtpi(0),
189 0 : fPtka(0),
190 0 : fETApic(0),
191 0 : fETAkac(0),
192 0 : fDecayer(0)
193 0 : {
194 : //
195 : // Standard constructor
196 : //
197 0 : fName="HIJINGpara";
198 0 : fTitle="HIJING Parametrisation Particle Generator";
199 0 : SetCutVertexZ();
200 0 : SetPtRange();
201 0 : }
202 :
203 : //_____________________________________________________________________________
204 0 : AliGenHIJINGpara::~AliGenHIJINGpara()
205 0 : {
206 : //
207 : // Standard destructor
208 : //
209 0 : delete fPtpi;
210 0 : delete fPtka;
211 0 : delete fETApic;
212 0 : delete fETAkac;
213 0 : }
214 :
215 : //_____________________________________________________________________________
216 : void AliGenHIJINGpara::Init()
217 : {
218 : //
219 : // Initialise the HIJING parametrisation
220 : //
221 0 : Float_t etaMin =-TMath::Log(TMath::Tan(
222 0 : TMath::Min((Double_t)fThetaMax/2,TMath::Pi()/2-1.e-10)));
223 0 : Float_t etaMax = -TMath::Log(TMath::Tan(
224 0 : TMath::Max((Double_t)fThetaMin/2,1.e-10)));
225 0 : fPtpi = new TF1("ptpi",&ptpi,0,20,0);
226 0 : gROOT->GetListOfFunctions()->Remove(fPtpi);
227 0 : fPtka = new TF1("ptka",&ptka,0,20,0);
228 0 : gROOT->GetListOfFunctions()->Remove(fPtka);
229 0 : fPtpi->SetNpx(1000);
230 0 : fPtka->SetNpx(1000);
231 0 : fETApic = new TF1("etapic",&etapic,etaMin,etaMax,0);
232 0 : gROOT->GetListOfFunctions()->Remove(fETApic);
233 0 : fETAkac = new TF1("etakac",&etakac,etaMin,etaMax,0);
234 0 : gROOT->GetListOfFunctions()->Remove(fETAkac);
235 :
236 0 : TF1 etaPic0("etaPic0",&etapic,-7,7,0);
237 0 : TF1 etaKac0("etaKac0",&etakac,-7,7,0);
238 :
239 0 : TF1 ptPic0("ptPic0",&ptpi,0.,15.,0);
240 0 : TF1 ptKac0("ptKac0",&ptka,0.,15.,0);
241 :
242 0 : Float_t intETApi = etaPic0.Integral(-0.5, 0.5);
243 0 : Float_t intETAka = etaKac0.Integral(-0.5, 0.5);
244 0 : Float_t scalePi = 7316/(intETApi/1.5);
245 0 : Float_t scaleKa = 684/(intETAka/2.0);
246 :
247 : // Fraction of events corresponding to the selected pt-range
248 0 : Float_t intPt = (0.877*ptPic0.Integral(0, 15)+
249 0 : 0.123*ptKac0.Integral(0, 15));
250 0 : Float_t intPtSel = (0.877*ptPic0.Integral(fPtMin, fPtMax)+
251 0 : 0.123*ptKac0.Integral(fPtMin, fPtMax));
252 0 : Float_t ptFrac = intPtSel/intPt;
253 :
254 : // Fraction of events corresponding to the selected eta-range
255 0 : Float_t intETASel = (scalePi*etaPic0.Integral(etaMin, etaMax)+
256 0 : scaleKa*etaKac0.Integral(etaMin, etaMax));
257 : // Fraction of events corresponding to the selected phi-range
258 0 : Float_t phiFrac = (fPhiMax-fPhiMin)/2/TMath::Pi();
259 :
260 :
261 0 : fParentWeight = (intETASel*ptFrac*phiFrac) / Float_t(fNpart);
262 :
263 0 : if (fAnalog != 0) {
264 0 : fPtWgtPi = (fPtMax - fPtMin) / fPtpi->Integral(0., 20.);
265 0 : fPtWgtKa = (fPtMax - fPtMin) / fPtka->Integral(0., 20.);
266 0 : fParentWeight = (intETASel*phiFrac) / Float_t(fNpart);
267 0 : }
268 :
269 :
270 0 : AliInfo(Form("The number of particles in the selected kinematic region corresponds to %f percent of a full event",
271 : 100./ fParentWeight));
272 :
273 : // Issue warning message if etaMin or etaMax are outside the alowed range
274 : // of the parametrization
275 0 : if (etaMin < -8.001 || etaMax > 8.001) {
276 0 : AliWarning("\nYOU ARE USING THE PARAMETERISATION OUTSIDE ");
277 0 : AliWarning("THE ALLOWED PSEUDORAPIDITY RANGE (-8. - 8.)");
278 0 : AliWarning(Form("YOUR LIMITS: %f %f \n ", etaMin, etaMax));
279 : }
280 : //
281 : //
282 0 : if (fPi0Decays && TVirtualMC::GetMC())
283 0 : fDecayer = TVirtualMC::GetMC()->GetDecayer();
284 :
285 0 : if (fPi0Decays)
286 : {
287 0 : fDecayer->SetForceDecay(kNeutralPion);
288 0 : fDecayer->Init();
289 : }
290 :
291 0 : }
292 :
293 :
294 : //_____________________________________________________________________________
295 : void AliGenHIJINGpara::Generate()
296 : {
297 : //
298 : // Generate one trigger
299 : //
300 :
301 :
302 : const Float_t kRaKpic=0.14;
303 : const Float_t kBorne=1/(1+kRaKpic);
304 0 : Float_t polar[3]= {0,0,0};
305 : //
306 : const Int_t kPions[3] = {kPi0, kPiPlus, kPiMinus};
307 : const Int_t kKaons[4] = {kK0Long, kK0Short, kKPlus, kKMinus};
308 : //
309 0 : Float_t origin[3];
310 : Float_t time;
311 : Float_t pt, pl, ptot, wgt;
312 : Float_t phi, theta;
313 0 : Float_t p[3];
314 : Int_t i, part, j;
315 : //
316 : TF1 *ptf;
317 : TF1 *etaf;
318 : //
319 0 : Float_t random[6];
320 : //
321 0 : for (j=0;j<3;j++) origin[j]=fOrigin[j];
322 0 : time = fTimeOrigin;
323 :
324 0 : if(fVertexSmear == kPerEvent) {
325 0 : Vertex();
326 0 : for (j=0; j < 3; j++) origin[j] = fVertex[j];
327 0 : time = fTime;
328 0 : } // if kPerEvent
329 0 : TArrayF eventVertex;
330 0 : eventVertex.Set(3);
331 0 : eventVertex[0] = origin[0];
332 0 : eventVertex[1] = origin[1];
333 0 : eventVertex[2] = origin[2];
334 : Float_t eventTime = time;
335 :
336 0 : for(i=0;i<fNpart;i++) {
337 : while(1) {
338 0 : Rndm(random,4);
339 0 : if(random[0]<kBorne) {
340 0 : part=kPions[Int_t (random[1]*3)];
341 0 : ptf=fPtpi;
342 0 : etaf=fETApic;
343 0 : wgt = fPtWgtPi;
344 0 : } else {
345 0 : part=kKaons[Int_t (random[1]*4)];
346 0 : ptf=fPtka;
347 0 : etaf=fETAkac;
348 0 : wgt = fPtWgtKa;
349 : }
350 0 : phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
351 0 : theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
352 0 : if(theta<fThetaMin || theta>fThetaMax) continue;
353 :
354 0 : if (fAnalog == 0) {
355 0 : pt = ptf->GetRandom();
356 0 : } else {
357 0 : pt = fPtMin + random[3] * (fPtMax - fPtMin);
358 : }
359 :
360 :
361 0 : pl=pt/TMath::Tan(theta);
362 0 : ptot=TMath::Sqrt(pt*pt+pl*pl);
363 0 : if(ptot<fPMin || ptot>fPMax) continue;
364 0 : p[0]=pt*TMath::Cos(phi);
365 0 : p[1]=pt*TMath::Sin(phi);
366 0 : p[2]=pl;
367 0 : if(fVertexSmear==kPerTrack) {
368 0 : Rndm(random,6);
369 0 : for (j=0;j<3;j++) {
370 0 : origin[j]=fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
371 0 : TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
372 : }
373 0 : Rndm(random,2);
374 0 : time = fTimeOrigin + fOsigma[2]/TMath::Ccgs()*
375 0 : TMath::Cos(2*random[0]*TMath::Pi())*
376 0 : TMath::Sqrt(-2*TMath::Log(random[1]));
377 0 : }
378 :
379 0 : if (fAnalog == 0) {
380 0 : wgt = fParentWeight;
381 0 : } else {
382 0 : wgt *= (fParentWeight * ptf->Eval(pt));
383 : }
384 :
385 :
386 0 : if (part == kPi0 && fPi0Decays){
387 : //
388 : // Decay pi0 if requested
389 0 : PushTrack(0,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);
390 0 : KeepTrack(fNt);
391 0 : DecayPi0(origin, p, time);
392 : } else {
393 : // printf("fNt %d", fNt);
394 0 : PushTrack(fTrackIt,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);
395 :
396 0 : KeepTrack(fNt);
397 : }
398 :
399 : break;
400 : }
401 0 : SetHighWaterMark(fNt);
402 : }
403 : //
404 :
405 : // Header
406 0 : AliGenEventHeader* header = new AliGenEventHeader("HIJINGparam");
407 : // Event Vertex
408 0 : header->SetPrimaryVertex(eventVertex);
409 0 : header->SetInteractionTime(eventTime);
410 0 : header->SetNProduced(fNpartProd);
411 0 : if (fContainer) {
412 0 : header->SetName(fName);
413 0 : fContainer->AddHeader(header);
414 : } else {
415 0 : gAlice->SetGenEventHeader(header);
416 : }
417 0 : }
418 :
419 : void AliGenHIJINGpara::SetPtRange(Float_t ptmin, Float_t ptmax) {
420 0 : AliGenerator::SetPtRange(ptmin, ptmax);
421 0 : }
422 :
423 : void AliGenHIJINGpara::DecayPi0(Float_t* orig, Float_t * p, Float_t time)
424 : {
425 : //
426 : // Decay the pi0
427 : // and put decay products on the stack
428 : //
429 : static TClonesArray *particles;
430 0 : if(!particles) particles = new TClonesArray("TParticle",1000);
431 : //
432 0 : const Float_t kMass = TDatabasePDG::Instance()->GetParticle(kPi0)->Mass();
433 0 : Float_t e = TMath::Sqrt(p[0] * p[0] + p[1] * p[1] + p[2] * p[2]+ kMass * kMass);
434 : //
435 : // Decay the pi0
436 0 : TLorentzVector pmom(p[0], p[1], p[2], e);
437 0 : fDecayer->Decay(kPi0, &pmom);
438 :
439 : //
440 : // Put decay particles on the stack
441 : //
442 0 : Float_t polar[3] = {0., 0., 0.};
443 0 : Int_t np = fDecayer->ImportParticles(particles);
444 0 : fNpartProd += (np-1);
445 0 : Int_t nt = 0;
446 0 : for (Int_t i = 1; i < np; i++)
447 : {
448 0 : TParticle* iParticle = (TParticle *) particles->At(i);
449 0 : p[0] = iParticle->Px();
450 0 : p[1] = iParticle->Py();
451 0 : p[2] = iParticle->Pz();
452 0 : Int_t part = iParticle->GetPdgCode();
453 :
454 0 : PushTrack(fTrackIt, fNt, part, p, orig, polar, time, kPDecay, nt, fParentWeight);
455 0 : KeepTrack(nt);
456 : }
457 0 : fNt = nt;
458 0 : }
459 :
460 : void AliGenHIJINGpara::Draw( const char * /*opt*/)
461 : {
462 : //
463 : // Draw the pT and y Distributions
464 : //
465 0 : TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
466 0 : c0->Divide(2,1);
467 0 : c0->cd(1);
468 0 : fPtpi->Draw();
469 0 : fPtpi->GetHistogram()->SetXTitle("p_{T} (GeV)");
470 0 : c0->cd(2);
471 0 : fPtka->Draw();
472 0 : fPtka->GetHistogram()->SetXTitle("p_{T} (GeV)");
473 :
474 0 : }
|