]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenParam.cxx
Cleaning up warnings (Sun)
[u/mrichter/AliRoot.git] / EVGEN / AliGenParam.cxx
CommitLineData
4c039060 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
88cb7938 16/* $Id$ */
6d4dd317 17
18// Class to generate particles from using paramtrized pT and y distributions.
19// Distributions are obtained from pointer to object of type AliGenLib.
20// (For example AliGenMUONlib)
21// Decays are performed using Pythia.
22// andreas.morsch@cern.ch
23
fe4da5cc 24#include "AliGenParam.h"
2904363f 25#include "AliDecayer.h"
fe4da5cc 26#include "AliGenMUONlib.h"
27#include "AliRun.h"
1578254f 28#include <TParticle.h>
fa480368 29#include <TParticlePDG.h>
30#include <TDatabasePDG.h>
31#include <TLorentzVector.h>
5d12ce38 32#include <TVirtualMC.h>
fa480368 33
f87cfe57 34#include <TF1.h>
5bd39445 35#include <TCanvas.h>
36#include <TH1.h>
5d12ce38 37#include "AliMC.h"
fe4da5cc 38
39ClassImp(AliGenParam)
40
41//------------------------------------------------------------
42
43 //Begin_Html
44 /*
1439f98e 45 <img src="picts/AliGenParam.gif">
fe4da5cc 46 */
47 //End_Html
48
49//____________________________________________________________
50//____________________________________________________________
51AliGenParam::AliGenParam()
fe4da5cc 52{
b22ee262 53// Deafault constructor
54 fPtPara = 0;
55 fYPara = 0;
34f60c01 56 fParam = 0;
57 fAnalog = kAnalog;
b22ee262 58 SetDeltaPt();
2685bf00 59 fDecayer = 0;
8b31bfac 60
61
b22ee262 62}
63
fff02fee 64AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library, Int_t param, char* tname):AliGenMC(npart)
b22ee262 65{
66// Constructor using number of particles parameterisation id and library
8b31bfac 67 fName = "Param";
68 fTitle= "Particle Generator using pT and y parameterisation";
b22ee262 69
70 fPtParaFunc = Library->GetPt(param, tname);
71 fYParaFunc = Library->GetY (param, tname);
72 fIpParaFunc = Library->GetIp(param, tname);
73
74 fPtPara = 0;
75 fYPara = 0;
76 fParam = param;
34f60c01 77 fAnalog = kAnalog;
b22ee262 78 SetForceDecay();
b22ee262 79 SetDeltaPt();
fe4da5cc 80}
81
82//____________________________________________________________
886b6f73 83
fff02fee 84AliGenParam::AliGenParam(Int_t npart, Int_t param, char* tname):AliGenMC(npart)
fe4da5cc 85{
b22ee262 86// Constructor using parameterisation id and number of particles
8b31bfac 87//
88 fName = "Param";
89 fTitle= "Particle Generator using pT and y parameterisation";
90
675e9664 91 AliGenLib* pLibrary = new AliGenMUONlib();
b22ee262 92
675e9664 93 fPtParaFunc = pLibrary->GetPt(param, tname);
94 fYParaFunc = pLibrary->GetY (param, tname);
95 fIpParaFunc = pLibrary->GetIp(param, tname);
b22ee262 96
97 fPtPara = 0;
98 fYPara = 0;
99 fParam = param;
34f60c01 100 fAnalog = kAnalog;
b22ee262 101 fChildSelect.Set(5);
102 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
103 SetForceDecay();
104 SetCutOnChild();
105 SetChildMomentumRange();
106 SetChildPtRange();
107 SetChildPhiRange();
108 SetChildThetaRange();
109 SetDeltaPt();
fe4da5cc 110}
111
34f60c01 112AliGenParam::AliGenParam(Int_t npart, Int_t param,
886b6f73 113 Double_t (*PtPara) (Double_t*, Double_t*),
114 Double_t (*YPara ) (Double_t* ,Double_t*),
65fb704d 115 Int_t (*IpPara) (TRandom *))
fff02fee 116 :AliGenMC(npart)
886b6f73 117{
749070b6 118// Constructor
28337bc1 119// Gines Martinez 1/10/99
8b31bfac 120 fName = "Param";
121 fTitle= "Particle Generator using pT and y parameterisation";
122
886b6f73 123 fPtParaFunc = PtPara;
124 fYParaFunc = YPara;
125 fIpParaFunc = IpPara;
28337bc1 126//
886b6f73 127 fPtPara = 0;
128 fYPara = 0;
129 fParam = param;
34f60c01 130 fAnalog = kAnalog;
886b6f73 131 fChildSelect.Set(5);
132 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
133 SetForceDecay();
134 SetCutOnChild();
749070b6 135 SetChildMomentumRange();
136 SetChildPtRange();
137 SetChildPhiRange();
138 SetChildThetaRange();
139 SetDeltaPt();
886b6f73 140}
141
f87cfe57 142
6d4dd317 143AliGenParam::AliGenParam(const AliGenParam & Param)
198bb1c7 144 :AliGenMC(Param)
f87cfe57 145{
198bb1c7 146// Copy constructor
6d4dd317 147 Param.Copy(*this);
f87cfe57 148}
149
fe4da5cc 150//____________________________________________________________
151AliGenParam::~AliGenParam()
152{
749070b6 153// Destructor
fe4da5cc 154 delete fPtPara;
155 delete fYPara;
156}
157
158//____________________________________________________________
159void AliGenParam::Init()
160{
749070b6 161// Initialisation
fa480368 162
2904363f 163 if (gMC) fDecayer = gMC->GetDecayer();
fe4da5cc 164 //Begin_Html
165 /*
1439f98e 166 <img src="picts/AliGenParam.gif">
fe4da5cc 167 */
168 //End_Html
169
749070b6 170 fPtPara = new TF1("Pt-Parametrization",fPtParaFunc,fPtMin,fPtMax,0);
171// Set representation precision to 10 MeV
172 Int_t npx= Int_t((fPtMax-fPtMin)/fDeltaPt);
173
174 fPtPara->SetNpx(npx);
175
2a336e15 176 fYPara = new TF1("Y-Parametrization",fYParaFunc,fYMin,fYMax,0);
177 TF1 ptPara("Pt-Parametrization(0,15)",fPtParaFunc,0,15,0);
178 TF1 yPara("Y-Parametrization(-6,6)",fYParaFunc,-6,6,0);
b7601ac4 179
fe4da5cc 180//
181// dN/dy| y=0
749070b6 182 Double_t y1=0;
183 Double_t y2=0;
184
185 fdNdy0=fYParaFunc(&y1,&y2);
fe4da5cc 186//
187// Integral over generation region
2a336e15 188 Float_t intYS = yPara.Integral(fYMin, fYMax);
189 Float_t intPt0 = ptPara.Integral(0,15);
190 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax);
749070b6 191 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
34f60c01 192 if (fAnalog == kAnalog) {
749070b6 193 fYWgt = intYS/fdNdy0;
194 fPtWgt = intPtS/intPt0;
195 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
196 } else {
197 fYWgt = intYS/fdNdy0;
198 fPtWgt = (fPtMax-fPtMin)/intPt0;
199 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
200 }
fe4da5cc 201//
202// particle decay related initialization
00d6ce7d 203 fDecayer->SetForceDecay(fForceDecay);
fa480368 204 fDecayer->Init();
00d6ce7d 205
fe4da5cc 206//
fff02fee 207 AliGenMC::Init();
fe4da5cc 208}
209
210//____________________________________________________________
211void AliGenParam::Generate()
212{
28337bc1 213//
214// Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
215// Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
216// antineutrons in the the desired theta, phi and momentum windows;
217// Gaussian smearing on the vertex is done if selected.
cc5d764c 218// The decay of heavy mesons is done using lujet,
219// and the childern particle are tracked by GEANT
220// However, light mesons are directly tracked by GEANT
221// setting fForceDecay = nodecay (SetForceDecay(nodecay))
28337bc1 222//
771bbe45 223//
224// Reinitialize decayer
225 fDecayer->Init();
226//
28337bc1 227 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
228 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
229 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
230 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
231 Float_t p[3], pc[3],
fff02fee 232 och[3]; // Momentum, polarisation and origin of the children particles from lujet
2d7a47be 233 Double_t ty, xmt;
fff02fee 234 Int_t nt, i, j;
fe4da5cc 235 Float_t wgtp, wgtch;
236 Double_t dummy;
09fd3ea2 237 static TClonesArray *particles;
fe4da5cc 238 //
fff02fee 239 if(!particles) particles = new TClonesArray("TParticle",1000);
8b31bfac 240
349be858 241 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
fe4da5cc 242 //
243 Float_t random[6];
28337bc1 244
245// Calculating vertex position per event
fe4da5cc 246 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
aee8290b 247 if(fVertexSmear==kPerEvent) {
d9ea0e3b 248 Vertex();
249 for (j=0;j<3;j++) origin0[j]=fVertex[j];
fe4da5cc 250 }
d9ea0e3b 251
21aaa175 252 Int_t ipa=0;
fff02fee 253
28337bc1 254// Generating fNpart particles
21aaa175 255 while (ipa<fNpart) {
cc5d764c 256 while(1) {
fe4da5cc 257//
258// particle type
65fb704d 259 Int_t iPart = fIpParaFunc(fRandom);
fa480368 260 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
675e9664 261 TParticlePDG *particle = pDataBase->GetParticle(iPart);
fa480368 262 Float_t am = particle->Mass();
8b31bfac 263
65fb704d 264 Rndm(random,2);
fe4da5cc 265//
266// phi
267 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
268//
269// y
2d7a47be 270 ty = TMath::TanH(fYPara->GetRandom());
fe4da5cc 271//
272// pT
34f60c01 273 if (fAnalog == kAnalog) {
fe4da5cc 274 pt=fPtPara->GetRandom();
275 wgtp=fParentWeight;
276 wgtch=fChildWeight;
277 } else {
278 pt=fPtMin+random[1]*(fPtMax-fPtMin);
279 Double_t ptd=pt;
280 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
281 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
282 }
283 xmt=sqrt(pt*pt+am*am);
2d7a47be 284 if (TMath::Abs(ty)==1.) {
285 ty=0.;
286 Fatal("AliGenParam",
287 "Division by 0: Please check you rapidity range !");
288 }
289
fe4da5cc 290 pl=xmt*ty/sqrt(1.-ty*ty);
291 theta=TMath::ATan2(pt,pl);
cc5d764c 292// Cut on theta
fe4da5cc 293 if(theta<fThetaMin || theta>fThetaMax) continue;
294 ptot=TMath::Sqrt(pt*pt+pl*pl);
cc5d764c 295// Cut on momentum
fe4da5cc 296 if(ptot<fPMin || ptot>fPMax) continue;
fff02fee 297//
fe4da5cc 298 p[0]=pt*TMath::Cos(phi);
299 p[1]=pt*TMath::Sin(phi);
300 p[2]=pl;
aee8290b 301 if(fVertexSmear==kPerTrack) {
65fb704d 302 Rndm(random,6);
fe4da5cc 303 for (j=0;j<3;j++) {
304 origin0[j]=
305 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
306 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
307 }
308 }
28337bc1 309
310// Looking at fForceDecay :
cc5d764c 311// if fForceDecay != none Primary particle decays using
312// AliPythia and children are tracked by GEANT
28337bc1 313//
cc5d764c 314// if fForceDecay == none Primary particle is tracked by GEANT
315// (In the latest, make sure that GEANT actually does all the decays you want)
fe4da5cc 316//
fff02fee 317
34f60c01 318 if (fForceDecay != kNoDecay) {
28337bc1 319// Using lujet to decay particle
886b6f73 320 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
fa480368 321 TLorentzVector pmom(p[0], p[1], p[2], energy);
322 fDecayer->Decay(iPart,&pmom);
fe4da5cc 323//
cc5d764c 324// select decay particles
fa480368 325 Int_t np=fDecayer->ImportParticles(particles);
1242532d 326
327 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
328 if (fGeometryAcceptance)
329 if (!CheckAcceptanceGeometry(np,particles)) continue;
886b6f73 330 Int_t ncsel=0;
fff02fee 331 Int_t* pFlag = new Int_t[np];
332 Int_t* pParent = new Int_t[np];
333 Int_t* pSelected = new Int_t[np];
334 Int_t* trackIt = new Int_t[np];
335
336 for (i=0; i<np; i++) {
337 pFlag[i] = 0;
338 pSelected[i] = 0;
339 pParent[i] = -1;
340 }
341
6ba00c52 342 if (np >1) {
343 TParticle* iparticle = (TParticle *) particles->At(0);
fff02fee 344 Int_t ipF, ipL;
345 for (i = 1; i<np ; i++) {
346 trackIt[i] = 1;
6ba00c52 347 iparticle = (TParticle *) particles->At(i);
348 Int_t kf = iparticle->GetPdgCode();
fff02fee 349 Int_t ks = iparticle->GetStatusCode();
350// flagged particle
351
352 if (pFlag[i] == 1) {
fff02fee 353 ipF = iparticle->GetFirstDaughter();
354 ipL = iparticle->GetLastDaughter();
355 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
356 continue;
357 }
358
359// flag decay products of particles with long life-time (c tau > .3 mum)
360
361 if (ks != 1) {
2d7a47be 362// TParticlePDG *particle = pDataBase->GetParticle(kf);
fff02fee 363
364 Double_t lifeTime = fDecayer->GetLifetime(kf);
2d7a47be 365// Double_t mass = particle->Mass();
366// Double_t width = particle->Width();
47fc6bd5 367 if (lifeTime > (Double_t) fMaxLifeTime) {
fff02fee 368 ipF = iparticle->GetFirstDaughter();
369 ipL = iparticle->GetLastDaughter();
370 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
371 } else{
372 trackIt[i] = 0;
373 pSelected[i] = 1;
374 }
375 } // ks==1 ?
fe4da5cc 376//
377// children
47fc6bd5 378
fff02fee 379 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
6ba00c52 380 {
6ba00c52 381 if (fCutOnChild) {
fff02fee 382 pc[0]=iparticle->Px();
383 pc[1]=iparticle->Py();
384 pc[2]=iparticle->Pz();
385 Bool_t childok = KinematicSelection(iparticle, 1);
386 if(childok) {
387 pSelected[i] = 1;
6ba00c52 388 ncsel++;
389 } else {
390 ncsel=-1;
391 break;
392 } // child kine cuts
393 } else {
fff02fee 394 pSelected[i] = 1;
886b6f73 395 ncsel++;
6ba00c52 396 } // if child selection
397 } // select muon
398 } // decay particle loop
399 } // if decay products
400
886b6f73 401 Int_t iparent;
402 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
403 ipa++;
21aaa175 404//
fff02fee 405// Parent
642f15cf 406 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
fff02fee 407 pParent[0] = nt;
a99cf51f 408 KeepTrack(nt);
fff02fee 409//
410// Decay Products
411//
412 for (i = 1; i < np; i++) {
413 if (pSelected[i]) {
414 TParticle* iparticle = (TParticle *) particles->At(i);
415 Int_t kf = iparticle->GetPdgCode();
416 Int_t ipa = iparticle->GetFirstMother()-1;
417
418 och[0] = origin0[0]+iparticle->Vx()/10;
419 och[1] = origin0[1]+iparticle->Vy()/10;
420 och[2] = origin0[2]+iparticle->Vz()/10;
421 pc[0] = iparticle->Px();
422 pc[1] = iparticle->Py();
423 pc[2] = iparticle->Pz();
424
425 if (ipa > -1) {
426 iparent = pParent[ipa];
427 } else {
428 iparent = -1;
429 }
1242532d 430
642f15cf 431 PushTrack(fTrackIt*trackIt[i], iparent, kf,
fff02fee 432 pc, och, polar,
433 0, kPDecay, nt, wgtch);
434 pParent[i] = nt;
a99cf51f 435 KeepTrack(nt);
8b31bfac 436 } // Selected
437 } // Particle loop
28337bc1 438 } // Decays by Lujet
8b31bfac 439 particles->Clear();
fff02fee 440 if (pFlag) delete[] pFlag;
441 if (pParent) delete[] pParent;
442 if (pSelected) delete[] pSelected;
8b31bfac 443 if (trackIt) delete[] trackIt;
21aaa175 444 } // kinematic selection
28337bc1 445 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
446 {
5d12ce38 447 gAlice->GetMCApp()->
642f15cf 448 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
28337bc1 449 ipa++;
450 }
fe4da5cc 451 break;
21aaa175 452 } // while
fe4da5cc 453 } // event loop
a99cf51f 454 SetHighWaterMark(nt);
fe4da5cc 455}
456
dc1d768c 457void AliGenParam::Draw( const char * /*opt*/)
5bd39445 458{
459 //
460 // Draw the pT and y Distributions
461 //
462 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
463 c0->Divide(2,1);
464 c0->cd(1);
465 fPtPara->Draw();
466 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
467 c0->cd(2);
468 fYPara->Draw();
469 fYPara->GetHistogram()->SetXTitle("y");
470}
471
f87cfe57 472AliGenParam& AliGenParam::operator=(const AliGenParam& rhs)
473{
474// Assignment operator
198bb1c7 475 rhs.Copy(*this);
f87cfe57 476 return *this;
477}
478
fe4da5cc 479
480