]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenParam.cxx
Use CalcMass()
[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
16/*
17$Log$
1242532d 18Revision 1.37 2002/04/17 10:20:44 morsch
19Coding Rule violations corrected.
20
6d4dd317 21Revision 1.36 2002/02/08 16:50:50 morsch
22Add name and title in constructor.
23
8b31bfac 24Revision 1.35 2002/01/21 10:02:40 morsch
25ty is Double_t
26Abort if too high rapidity causes numerical paroblem. User has to specify
27meaningful y-range.
28
2d7a47be 29Revision 1.34 2001/11/27 13:13:07 morsch
30Maximum lifetime for long-lived particles to be put on the stack is parameter.
31It can be set via SetMaximumLifetime(..).
32
47fc6bd5 33Revision 1.33 2001/10/21 18:35:56 hristov
34Several pointers were set to zero in the default constructors to avoid memory management problems
35
2685bf00 36Revision 1.32 2001/07/27 17:09:36 morsch
37Use local SetTrack, KeepTrack and SetHighWaterMark methods
38to delegate either to local stack or to stack owned by AliRun.
39(Piotr Skowronski, A.M.)
40
a99cf51f 41Revision 1.31 2001/07/13 10:58:54 morsch
42- Some coded moved to AliGenMC
43- Improved handling of secondary vertices.
44
fff02fee 45Revision 1.30 2001/06/15 07:55:04 morsch
46Put only first generation decay products on the stack.
47
6ba00c52 48Revision 1.29 2001/03/27 10:58:41 morsch
49Initialize decayer before generation. Important if run inside cocktail.
50
771bbe45 51Revision 1.28 2001/03/09 13:01:41 morsch
52- enum constants for paramterisation type (particle family) moved to AliGen*lib.h
53- use AliGenGSIlib::kUpsilon, AliGenPHOSlib::kEtaPrime to access the constants
54
34f60c01 55Revision 1.27 2001/02/02 15:21:10 morsch
56Set high water mark after last particle.
57Use Vertex() method for Vertex.
58
d9ea0e3b 59Revision 1.26 2000/12/21 16:24:06 morsch
60Coding convention clean-up
61
675e9664 62Revision 1.25 2000/11/30 07:12:50 alibrary
63Introducing new Rndm and QA classes
64
65fb704d 65Revision 1.24 2000/10/18 19:11:27 hristov
66Division by zero fixed
67
919c2096 68Revision 1.23 2000/10/02 21:28:06 fca
69Removal of useless dependecies via forward declarations
70
94de3818 71Revision 1.22 2000/09/12 14:14:55 morsch
72Call fDecayer->ForceDecay() at the beginning of Generate().
73
00d6ce7d 74Revision 1.21 2000/09/08 15:39:01 morsch
75Handle the case fForceDecay=all during the generation, i.e. select all secondaries.
76
5ab0acc9 77Revision 1.20 2000/09/06 14:35:44 morsch
78Use AliDecayerPythia for particle decays.
79
fa480368 80Revision 1.19 2000/07/11 18:24:56 fca
81Coding convention corrections + few minor bug fixes
82
aee8290b 83Revision 1.18 2000/06/29 21:08:27 morsch
84All paramatrisation libraries derive from the pure virtual base class AliGenLib.
85This allows to pass a pointer to a library directly to AliGenParam and avoids the
86use of function pointers in Config.C.
87
b22ee262 88Revision 1.17 2000/06/09 20:33:30 morsch
89All coding rule violations except RS3 corrected
90
f87cfe57 91Revision 1.16 2000/05/02 07:51:31 morsch
92- Control precision of pT sampling TF1::SetNpx(..)
93- Correct initialisation of child-cuts in all constructors.
94- Most coding rule violations corrected.
95
749070b6 96Revision 1.15 2000/04/03 15:42:12 morsch
97Cuts on primary particles are separated from those on the decay products. Methods
98SetChildMomentumRange, SetChildPtRange, SetChildPhiRange, SetChildThetaRange added.
99
cc5d764c 100Revision 1.14 1999/11/09 07:38:48 fca
101Changes for compatibility with version 2.23 of ROOT
102
084c1b4a 103Revision 1.13 1999/11/04 11:30:31 fca
104Correct the logics for SetForceDecay
105
28337bc1 106Revision 1.12 1999/11/03 17:43:20 fca
107New version from G.Martinez & A.Morsch
108
886b6f73 109Revision 1.11 1999/09/29 09:24:14 fca
110Introduction of the Copyright and cvs Log
111
4c039060 112*/
113
6d4dd317 114
115
116// Class to generate particles from using paramtrized pT and y distributions.
117// Distributions are obtained from pointer to object of type AliGenLib.
118// (For example AliGenMUONlib)
119// Decays are performed using Pythia.
120// andreas.morsch@cern.ch
121
fe4da5cc 122#include "AliGenParam.h"
fa480368 123#include "AliDecayerPythia.h"
fe4da5cc 124#include "AliGenMUONlib.h"
125#include "AliRun.h"
1578254f 126#include <TParticle.h>
fa480368 127#include <TParticlePDG.h>
128#include <TDatabasePDG.h>
129#include <TLorentzVector.h>
130
f87cfe57 131#include <TF1.h>
fe4da5cc 132
133ClassImp(AliGenParam)
134
135//------------------------------------------------------------
136
137 //Begin_Html
138 /*
1439f98e 139 <img src="picts/AliGenParam.gif">
fe4da5cc 140 */
141 //End_Html
142
143//____________________________________________________________
144//____________________________________________________________
145AliGenParam::AliGenParam()
fe4da5cc 146{
b22ee262 147// Deafault constructor
148 fPtPara = 0;
149 fYPara = 0;
34f60c01 150 fParam = 0;
151 fAnalog = kAnalog;
b22ee262 152 SetDeltaPt();
675e9664 153//
154// Set random number generator
155 sRandom = fRandom;
2685bf00 156 fDecayer = 0;
8b31bfac 157
158
b22ee262 159}
160
fff02fee 161AliGenParam::AliGenParam(Int_t npart, AliGenLib * Library, Int_t param, char* tname):AliGenMC(npart)
b22ee262 162{
163// Constructor using number of particles parameterisation id and library
8b31bfac 164 fName = "Param";
165 fTitle= "Particle Generator using pT and y parameterisation";
b22ee262 166
167 fPtParaFunc = Library->GetPt(param, tname);
168 fYParaFunc = Library->GetY (param, tname);
169 fIpParaFunc = Library->GetIp(param, tname);
170
171 fPtPara = 0;
172 fYPara = 0;
173 fParam = param;
34f60c01 174 fAnalog = kAnalog;
b22ee262 175 SetForceDecay();
b22ee262 176 SetDeltaPt();
675e9664 177//
178// Set random number generator
179 sRandom = fRandom;
fe4da5cc 180}
181
182//____________________________________________________________
886b6f73 183
fff02fee 184AliGenParam::AliGenParam(Int_t npart, Int_t param, char* tname):AliGenMC(npart)
fe4da5cc 185{
b22ee262 186// Constructor using parameterisation id and number of particles
8b31bfac 187//
188 fName = "Param";
189 fTitle= "Particle Generator using pT and y parameterisation";
190
675e9664 191 AliGenLib* pLibrary = new AliGenMUONlib();
b22ee262 192
675e9664 193 fPtParaFunc = pLibrary->GetPt(param, tname);
194 fYParaFunc = pLibrary->GetY (param, tname);
195 fIpParaFunc = pLibrary->GetIp(param, tname);
b22ee262 196
197 fPtPara = 0;
198 fYPara = 0;
199 fParam = param;
34f60c01 200 fAnalog = kAnalog;
b22ee262 201 fChildSelect.Set(5);
202 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
203 SetForceDecay();
204 SetCutOnChild();
205 SetChildMomentumRange();
206 SetChildPtRange();
207 SetChildPhiRange();
208 SetChildThetaRange();
209 SetDeltaPt();
fe4da5cc 210}
211
34f60c01 212AliGenParam::AliGenParam(Int_t npart, Int_t param,
886b6f73 213 Double_t (*PtPara) (Double_t*, Double_t*),
214 Double_t (*YPara ) (Double_t* ,Double_t*),
65fb704d 215 Int_t (*IpPara) (TRandom *))
fff02fee 216 :AliGenMC(npart)
886b6f73 217{
749070b6 218// Constructor
28337bc1 219// Gines Martinez 1/10/99
8b31bfac 220 fName = "Param";
221 fTitle= "Particle Generator using pT and y parameterisation";
222
886b6f73 223 fPtParaFunc = PtPara;
224 fYParaFunc = YPara;
225 fIpParaFunc = IpPara;
28337bc1 226//
886b6f73 227 fPtPara = 0;
228 fYPara = 0;
229 fParam = param;
34f60c01 230 fAnalog = kAnalog;
886b6f73 231 fChildSelect.Set(5);
232 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
233 SetForceDecay();
234 SetCutOnChild();
749070b6 235 SetChildMomentumRange();
236 SetChildPtRange();
237 SetChildPhiRange();
238 SetChildThetaRange();
239 SetDeltaPt();
886b6f73 240}
241
f87cfe57 242
6d4dd317 243AliGenParam::AliGenParam(const AliGenParam & Param)
f87cfe57 244{
245// copy constructor
6d4dd317 246 Param.Copy(*this);
f87cfe57 247}
248
fe4da5cc 249//____________________________________________________________
250AliGenParam::~AliGenParam()
251{
749070b6 252// Destructor
fe4da5cc 253 delete fPtPara;
254 delete fYPara;
255}
256
257//____________________________________________________________
258void AliGenParam::Init()
259{
749070b6 260// Initialisation
fa480368 261
262 fDecayer = new AliDecayerPythia();
fe4da5cc 263 //Begin_Html
264 /*
1439f98e 265 <img src="picts/AliGenParam.gif">
fe4da5cc 266 */
267 //End_Html
268
749070b6 269 fPtPara = new TF1("Pt-Parametrization",fPtParaFunc,fPtMin,fPtMax,0);
270// Set representation precision to 10 MeV
271 Int_t npx= Int_t((fPtMax-fPtMin)/fDeltaPt);
272
273 fPtPara->SetNpx(npx);
274
275 fYPara = new TF1("Y -Parametrization",fYParaFunc,fYMin,fYMax,0);
276 TF1* ptPara = new TF1("Pt-Parametrization",fPtParaFunc,0,15,0);
277 TF1* yPara = new TF1("Y -Parametrization",fYParaFunc,-6,6,0);
b7601ac4 278
fe4da5cc 279//
280// dN/dy| y=0
749070b6 281 Double_t y1=0;
282 Double_t y2=0;
283
284 fdNdy0=fYParaFunc(&y1,&y2);
fe4da5cc 285//
286// Integral over generation region
749070b6 287 Float_t intYS = yPara ->Integral(fYMin, fYMax);
288 Float_t intPt0 = ptPara->Integral(0,15);
289 Float_t intPtS = ptPara->Integral(fPtMin,fPtMax);
290 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
34f60c01 291 if (fAnalog == kAnalog) {
749070b6 292 fYWgt = intYS/fdNdy0;
293 fPtWgt = intPtS/intPt0;
294 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
295 } else {
296 fYWgt = intYS/fdNdy0;
297 fPtWgt = (fPtMax-fPtMin)/intPt0;
298 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
299 }
fe4da5cc 300//
301// particle decay related initialization
00d6ce7d 302 fDecayer->SetForceDecay(fForceDecay);
fa480368 303 fDecayer->Init();
00d6ce7d 304
fe4da5cc 305//
fff02fee 306 AliGenMC::Init();
fe4da5cc 307}
308
309//____________________________________________________________
310void AliGenParam::Generate()
311{
28337bc1 312//
313// Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
314// Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
315// antineutrons in the the desired theta, phi and momentum windows;
316// Gaussian smearing on the vertex is done if selected.
cc5d764c 317// The decay of heavy mesons is done using lujet,
318// and the childern particle are tracked by GEANT
319// However, light mesons are directly tracked by GEANT
320// setting fForceDecay = nodecay (SetForceDecay(nodecay))
28337bc1 321//
771bbe45 322//
323// Reinitialize decayer
324 fDecayer->Init();
325//
28337bc1 326 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
327 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
328 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
329 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
330 Float_t p[3], pc[3],
fff02fee 331 och[3]; // Momentum, polarisation and origin of the children particles from lujet
2d7a47be 332 Double_t ty, xmt;
fff02fee 333 Int_t nt, i, j;
fe4da5cc 334 Float_t wgtp, wgtch;
335 Double_t dummy;
09fd3ea2 336 static TClonesArray *particles;
fe4da5cc 337 //
fff02fee 338 if(!particles) particles = new TClonesArray("TParticle",1000);
8b31bfac 339
675e9664 340 static TDatabasePDG *pDataBase = new TDatabasePDG();
341 if(!pDataBase) pDataBase = new TDatabasePDG();
fe4da5cc 342 //
343 Float_t random[6];
28337bc1 344
345// Calculating vertex position per event
fe4da5cc 346 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
aee8290b 347 if(fVertexSmear==kPerEvent) {
d9ea0e3b 348 Vertex();
349 for (j=0;j<3;j++) origin0[j]=fVertex[j];
fe4da5cc 350 }
d9ea0e3b 351
21aaa175 352 Int_t ipa=0;
fff02fee 353
28337bc1 354// Generating fNpart particles
21aaa175 355 while (ipa<fNpart) {
cc5d764c 356 while(1) {
fe4da5cc 357//
358// particle type
65fb704d 359 Int_t iPart = fIpParaFunc(fRandom);
fa480368 360 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
675e9664 361 TParticlePDG *particle = pDataBase->GetParticle(iPart);
fa480368 362 Float_t am = particle->Mass();
8b31bfac 363
65fb704d 364 Rndm(random,2);
fe4da5cc 365//
366// phi
367 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
368//
369// y
2d7a47be 370 ty = TMath::TanH(fYPara->GetRandom());
fe4da5cc 371//
372// pT
34f60c01 373 if (fAnalog == kAnalog) {
fe4da5cc 374 pt=fPtPara->GetRandom();
375 wgtp=fParentWeight;
376 wgtch=fChildWeight;
377 } else {
378 pt=fPtMin+random[1]*(fPtMax-fPtMin);
379 Double_t ptd=pt;
380 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
381 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
382 }
383 xmt=sqrt(pt*pt+am*am);
2d7a47be 384 if (TMath::Abs(ty)==1.) {
385 ty=0.;
386 Fatal("AliGenParam",
387 "Division by 0: Please check you rapidity range !");
388 }
389
fe4da5cc 390 pl=xmt*ty/sqrt(1.-ty*ty);
391 theta=TMath::ATan2(pt,pl);
cc5d764c 392// Cut on theta
fe4da5cc 393 if(theta<fThetaMin || theta>fThetaMax) continue;
394 ptot=TMath::Sqrt(pt*pt+pl*pl);
cc5d764c 395// Cut on momentum
fe4da5cc 396 if(ptot<fPMin || ptot>fPMax) continue;
fff02fee 397//
fe4da5cc 398 p[0]=pt*TMath::Cos(phi);
399 p[1]=pt*TMath::Sin(phi);
400 p[2]=pl;
aee8290b 401 if(fVertexSmear==kPerTrack) {
65fb704d 402 Rndm(random,6);
fe4da5cc 403 for (j=0;j<3;j++) {
404 origin0[j]=
405 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
406 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
407 }
408 }
28337bc1 409
410// Looking at fForceDecay :
cc5d764c 411// if fForceDecay != none Primary particle decays using
412// AliPythia and children are tracked by GEANT
28337bc1 413//
cc5d764c 414// if fForceDecay == none Primary particle is tracked by GEANT
415// (In the latest, make sure that GEANT actually does all the decays you want)
fe4da5cc 416//
fff02fee 417
34f60c01 418 if (fForceDecay != kNoDecay) {
28337bc1 419// Using lujet to decay particle
886b6f73 420 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
fa480368 421 TLorentzVector pmom(p[0], p[1], p[2], energy);
422 fDecayer->Decay(iPart,&pmom);
fe4da5cc 423//
cc5d764c 424// select decay particles
fa480368 425 Int_t np=fDecayer->ImportParticles(particles);
1242532d 426
427 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
428 if (fGeometryAcceptance)
429 if (!CheckAcceptanceGeometry(np,particles)) continue;
430
886b6f73 431 Int_t ncsel=0;
fff02fee 432 Int_t* pFlag = new Int_t[np];
433 Int_t* pParent = new Int_t[np];
434 Int_t* pSelected = new Int_t[np];
435 Int_t* trackIt = new Int_t[np];
436
437 for (i=0; i<np; i++) {
438 pFlag[i] = 0;
439 pSelected[i] = 0;
440 pParent[i] = -1;
441 }
442
6ba00c52 443 if (np >1) {
444 TParticle* iparticle = (TParticle *) particles->At(0);
fff02fee 445 Int_t ipF, ipL;
446 for (i = 1; i<np ; i++) {
447 trackIt[i] = 1;
6ba00c52 448 iparticle = (TParticle *) particles->At(i);
449 Int_t kf = iparticle->GetPdgCode();
fff02fee 450 Int_t ks = iparticle->GetStatusCode();
451// flagged particle
452
453 if (pFlag[i] == 1) {
fff02fee 454 ipF = iparticle->GetFirstDaughter();
455 ipL = iparticle->GetLastDaughter();
456 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
457 continue;
458 }
459
460// flag decay products of particles with long life-time (c tau > .3 mum)
461
462 if (ks != 1) {
2d7a47be 463// TParticlePDG *particle = pDataBase->GetParticle(kf);
fff02fee 464
465 Double_t lifeTime = fDecayer->GetLifetime(kf);
2d7a47be 466// Double_t mass = particle->Mass();
467// Double_t width = particle->Width();
47fc6bd5 468 if (lifeTime > (Double_t) fMaxLifeTime) {
fff02fee 469 ipF = iparticle->GetFirstDaughter();
470 ipL = iparticle->GetLastDaughter();
471 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
472 } else{
473 trackIt[i] = 0;
474 pSelected[i] = 1;
475 }
476 } // ks==1 ?
fe4da5cc 477//
478// children
47fc6bd5 479
fff02fee 480 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
6ba00c52 481 {
6ba00c52 482 if (fCutOnChild) {
fff02fee 483 pc[0]=iparticle->Px();
484 pc[1]=iparticle->Py();
485 pc[2]=iparticle->Pz();
486 Bool_t childok = KinematicSelection(iparticle, 1);
487 if(childok) {
488 pSelected[i] = 1;
6ba00c52 489 ncsel++;
490 } else {
491 ncsel=-1;
492 break;
493 } // child kine cuts
494 } else {
fff02fee 495 pSelected[i] = 1;
886b6f73 496 ncsel++;
6ba00c52 497 } // if child selection
498 } // select muon
499 } // decay particle loop
500 } // if decay products
501
886b6f73 502 Int_t iparent;
503 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
504 ipa++;
21aaa175 505//
fff02fee 506// Parent
a99cf51f 507 SetTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
fff02fee 508 pParent[0] = nt;
a99cf51f 509 KeepTrack(nt);
fff02fee 510//
511// Decay Products
512//
513 for (i = 1; i < np; i++) {
514 if (pSelected[i]) {
515 TParticle* iparticle = (TParticle *) particles->At(i);
516 Int_t kf = iparticle->GetPdgCode();
517 Int_t ipa = iparticle->GetFirstMother()-1;
518
519 och[0] = origin0[0]+iparticle->Vx()/10;
520 och[1] = origin0[1]+iparticle->Vy()/10;
521 och[2] = origin0[2]+iparticle->Vz()/10;
522 pc[0] = iparticle->Px();
523 pc[1] = iparticle->Py();
524 pc[2] = iparticle->Pz();
525
526 if (ipa > -1) {
527 iparent = pParent[ipa];
528 } else {
529 iparent = -1;
530 }
1242532d 531
a99cf51f 532 SetTrack(fTrackIt*trackIt[i], iparent, kf,
fff02fee 533 pc, och, polar,
534 0, kPDecay, nt, wgtch);
535 pParent[i] = nt;
a99cf51f 536 KeepTrack(nt);
8b31bfac 537 } // Selected
538 } // Particle loop
28337bc1 539 } // Decays by Lujet
8b31bfac 540 particles->Clear();
fff02fee 541 if (pFlag) delete[] pFlag;
542 if (pParent) delete[] pParent;
543 if (pSelected) delete[] pSelected;
8b31bfac 544 if (trackIt) delete[] trackIt;
21aaa175 545 } // kinematic selection
28337bc1 546 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
547 {
548 gAlice->
65fb704d 549 SetTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
28337bc1 550 ipa++;
551 }
fe4da5cc 552 break;
21aaa175 553 } // while
fe4da5cc 554 } // event loop
a99cf51f 555 SetHighWaterMark(nt);
fe4da5cc 556}
557
f87cfe57 558AliGenParam& AliGenParam::operator=(const AliGenParam& rhs)
559{
560// Assignment operator
561 return *this;
562}
563
fe4da5cc 564
565