]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenParam.cxx
Unique names for functions enforced.
[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
c5c3e4b0 169 char name[256];
170 sprintf(name, "pt-parameterisation for %s", GetName());
171
172 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
749070b6 173// Set representation precision to 10 MeV
c5c3e4b0 174 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
749070b6 175
176 fPtPara->SetNpx(npx);
c5c3e4b0 177
178 sprintf(name, "y-parameterisation for %s", GetName());
179 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
749070b6 180
c5c3e4b0 181 sprintf(name, "pt-for-%s", GetName());
182 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
183 sprintf(name, "y-for-%s", GetName());
184 TF1 yPara(name, fYParaFunc, -6, 6, 0);
b7601ac4 185
fe4da5cc 186//
187// dN/dy| y=0
749070b6 188 Double_t y1=0;
189 Double_t y2=0;
190
191 fdNdy0=fYParaFunc(&y1,&y2);
fe4da5cc 192//
193// Integral over generation region
2a336e15 194 Float_t intYS = yPara.Integral(fYMin, fYMax);
195 Float_t intPt0 = ptPara.Integral(0,15);
196 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax);
749070b6 197 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi();
34f60c01 198 if (fAnalog == kAnalog) {
749070b6 199 fYWgt = intYS/fdNdy0;
200 fPtWgt = intPtS/intPt0;
201 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
202 } else {
203 fYWgt = intYS/fdNdy0;
204 fPtWgt = (fPtMax-fPtMin)/intPt0;
205 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
206 }
fe4da5cc 207//
208// particle decay related initialization
00d6ce7d 209 fDecayer->SetForceDecay(fForceDecay);
fa480368 210 fDecayer->Init();
00d6ce7d 211
fe4da5cc 212//
fff02fee 213 AliGenMC::Init();
fe4da5cc 214}
215
216//____________________________________________________________
217void AliGenParam::Generate()
218{
28337bc1 219//
220// Generate 'npart' of light and heavy mesons (J/Psi, upsilon or phi, Pion,
221// Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
222// antineutrons in the the desired theta, phi and momentum windows;
223// Gaussian smearing on the vertex is done if selected.
cc5d764c 224// The decay of heavy mesons is done using lujet,
225// and the childern particle are tracked by GEANT
226// However, light mesons are directly tracked by GEANT
227// setting fForceDecay = nodecay (SetForceDecay(nodecay))
28337bc1 228//
771bbe45 229//
230// Reinitialize decayer
231 fDecayer->Init();
232//
28337bc1 233 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
234 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
235 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
236 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
237 Float_t p[3], pc[3],
fff02fee 238 och[3]; // Momentum, polarisation and origin of the children particles from lujet
2d7a47be 239 Double_t ty, xmt;
fff02fee 240 Int_t nt, i, j;
fe4da5cc 241 Float_t wgtp, wgtch;
242 Double_t dummy;
09fd3ea2 243 static TClonesArray *particles;
fe4da5cc 244 //
fff02fee 245 if(!particles) particles = new TClonesArray("TParticle",1000);
8b31bfac 246
349be858 247 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
fe4da5cc 248 //
249 Float_t random[6];
28337bc1 250
251// Calculating vertex position per event
fe4da5cc 252 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
aee8290b 253 if(fVertexSmear==kPerEvent) {
d9ea0e3b 254 Vertex();
255 for (j=0;j<3;j++) origin0[j]=fVertex[j];
fe4da5cc 256 }
d9ea0e3b 257
21aaa175 258 Int_t ipa=0;
fff02fee 259
28337bc1 260// Generating fNpart particles
21aaa175 261 while (ipa<fNpart) {
cc5d764c 262 while(1) {
fe4da5cc 263//
264// particle type
65fb704d 265 Int_t iPart = fIpParaFunc(fRandom);
fa480368 266 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
675e9664 267 TParticlePDG *particle = pDataBase->GetParticle(iPart);
fa480368 268 Float_t am = particle->Mass();
8b31bfac 269
65fb704d 270 Rndm(random,2);
fe4da5cc 271//
272// phi
273 phi=fPhiMin+random[0]*(fPhiMax-fPhiMin);
274//
275// y
2d7a47be 276 ty = TMath::TanH(fYPara->GetRandom());
fe4da5cc 277//
278// pT
34f60c01 279 if (fAnalog == kAnalog) {
fe4da5cc 280 pt=fPtPara->GetRandom();
281 wgtp=fParentWeight;
282 wgtch=fChildWeight;
283 } else {
284 pt=fPtMin+random[1]*(fPtMax-fPtMin);
285 Double_t ptd=pt;
286 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
287 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
288 }
289 xmt=sqrt(pt*pt+am*am);
2d7a47be 290 if (TMath::Abs(ty)==1.) {
291 ty=0.;
292 Fatal("AliGenParam",
293 "Division by 0: Please check you rapidity range !");
294 }
295
fe4da5cc 296 pl=xmt*ty/sqrt(1.-ty*ty);
297 theta=TMath::ATan2(pt,pl);
cc5d764c 298// Cut on theta
fe4da5cc 299 if(theta<fThetaMin || theta>fThetaMax) continue;
300 ptot=TMath::Sqrt(pt*pt+pl*pl);
cc5d764c 301// Cut on momentum
fe4da5cc 302 if(ptot<fPMin || ptot>fPMax) continue;
fff02fee 303//
fe4da5cc 304 p[0]=pt*TMath::Cos(phi);
305 p[1]=pt*TMath::Sin(phi);
306 p[2]=pl;
aee8290b 307 if(fVertexSmear==kPerTrack) {
65fb704d 308 Rndm(random,6);
fe4da5cc 309 for (j=0;j<3;j++) {
310 origin0[j]=
311 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
312 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
313 }
314 }
28337bc1 315
316// Looking at fForceDecay :
cc5d764c 317// if fForceDecay != none Primary particle decays using
318// AliPythia and children are tracked by GEANT
28337bc1 319//
cc5d764c 320// if fForceDecay == none Primary particle is tracked by GEANT
321// (In the latest, make sure that GEANT actually does all the decays you want)
fe4da5cc 322//
fff02fee 323
34f60c01 324 if (fForceDecay != kNoDecay) {
28337bc1 325// Using lujet to decay particle
886b6f73 326 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
fa480368 327 TLorentzVector pmom(p[0], p[1], p[2], energy);
328 fDecayer->Decay(iPart,&pmom);
fe4da5cc 329//
cc5d764c 330// select decay particles
fa480368 331 Int_t np=fDecayer->ImportParticles(particles);
1242532d 332
333 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
334 if (fGeometryAcceptance)
335 if (!CheckAcceptanceGeometry(np,particles)) continue;
886b6f73 336 Int_t ncsel=0;
fff02fee 337 Int_t* pFlag = new Int_t[np];
338 Int_t* pParent = new Int_t[np];
339 Int_t* pSelected = new Int_t[np];
340 Int_t* trackIt = new Int_t[np];
341
342 for (i=0; i<np; i++) {
343 pFlag[i] = 0;
344 pSelected[i] = 0;
345 pParent[i] = -1;
346 }
347
6ba00c52 348 if (np >1) {
349 TParticle* iparticle = (TParticle *) particles->At(0);
fff02fee 350 Int_t ipF, ipL;
351 for (i = 1; i<np ; i++) {
352 trackIt[i] = 1;
6ba00c52 353 iparticle = (TParticle *) particles->At(i);
354 Int_t kf = iparticle->GetPdgCode();
fff02fee 355 Int_t ks = iparticle->GetStatusCode();
356// flagged particle
357
358 if (pFlag[i] == 1) {
fff02fee 359 ipF = iparticle->GetFirstDaughter();
360 ipL = iparticle->GetLastDaughter();
361 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
362 continue;
363 }
364
365// flag decay products of particles with long life-time (c tau > .3 mum)
366
367 if (ks != 1) {
2d7a47be 368// TParticlePDG *particle = pDataBase->GetParticle(kf);
fff02fee 369
370 Double_t lifeTime = fDecayer->GetLifetime(kf);
2d7a47be 371// Double_t mass = particle->Mass();
372// Double_t width = particle->Width();
47fc6bd5 373 if (lifeTime > (Double_t) fMaxLifeTime) {
fff02fee 374 ipF = iparticle->GetFirstDaughter();
375 ipL = iparticle->GetLastDaughter();
376 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
377 } else{
378 trackIt[i] = 0;
379 pSelected[i] = 1;
380 }
381 } // ks==1 ?
fe4da5cc 382//
383// children
47fc6bd5 384
fff02fee 385 if (ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll && trackIt[i])
6ba00c52 386 {
6ba00c52 387 if (fCutOnChild) {
fff02fee 388 pc[0]=iparticle->Px();
389 pc[1]=iparticle->Py();
390 pc[2]=iparticle->Pz();
391 Bool_t childok = KinematicSelection(iparticle, 1);
392 if(childok) {
393 pSelected[i] = 1;
6ba00c52 394 ncsel++;
395 } else {
396 ncsel=-1;
397 break;
398 } // child kine cuts
399 } else {
fff02fee 400 pSelected[i] = 1;
886b6f73 401 ncsel++;
6ba00c52 402 } // if child selection
403 } // select muon
404 } // decay particle loop
405 } // if decay products
406
886b6f73 407 Int_t iparent;
408 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
409 ipa++;
21aaa175 410//
fff02fee 411// Parent
642f15cf 412 PushTrack(0, -1, iPart, p, origin0, polar, 0, kPPrimary, nt, wgtp);
fff02fee 413 pParent[0] = nt;
a99cf51f 414 KeepTrack(nt);
fff02fee 415//
416// Decay Products
417//
418 for (i = 1; i < np; i++) {
419 if (pSelected[i]) {
420 TParticle* iparticle = (TParticle *) particles->At(i);
421 Int_t kf = iparticle->GetPdgCode();
422 Int_t ipa = iparticle->GetFirstMother()-1;
423
424 och[0] = origin0[0]+iparticle->Vx()/10;
425 och[1] = origin0[1]+iparticle->Vy()/10;
426 och[2] = origin0[2]+iparticle->Vz()/10;
427 pc[0] = iparticle->Px();
428 pc[1] = iparticle->Py();
429 pc[2] = iparticle->Pz();
430
431 if (ipa > -1) {
432 iparent = pParent[ipa];
433 } else {
434 iparent = -1;
435 }
1242532d 436
642f15cf 437 PushTrack(fTrackIt*trackIt[i], iparent, kf,
fff02fee 438 pc, och, polar,
439 0, kPDecay, nt, wgtch);
440 pParent[i] = nt;
a99cf51f 441 KeepTrack(nt);
8b31bfac 442 } // Selected
443 } // Particle loop
28337bc1 444 } // Decays by Lujet
8b31bfac 445 particles->Clear();
fff02fee 446 if (pFlag) delete[] pFlag;
447 if (pParent) delete[] pParent;
448 if (pSelected) delete[] pSelected;
8b31bfac 449 if (trackIt) delete[] trackIt;
21aaa175 450 } // kinematic selection
28337bc1 451 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
452 {
5d12ce38 453 gAlice->GetMCApp()->
642f15cf 454 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,0,kPPrimary,nt,wgtp);
28337bc1 455 ipa++;
456 }
fe4da5cc 457 break;
21aaa175 458 } // while
fe4da5cc 459 } // event loop
a99cf51f 460 SetHighWaterMark(nt);
fe4da5cc 461}
462
dc1d768c 463void AliGenParam::Draw( const char * /*opt*/)
5bd39445 464{
465 //
466 // Draw the pT and y Distributions
467 //
468 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
469 c0->Divide(2,1);
470 c0->cd(1);
471 fPtPara->Draw();
472 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
473 c0->cd(2);
474 fYPara->Draw();
475 fYPara->GetHistogram()->SetXTitle("y");
476}
477
f87cfe57 478AliGenParam& AliGenParam::operator=(const AliGenParam& rhs)
479{
480// Assignment operator
198bb1c7 481 rhs.Copy(*this);
f87cfe57 482 return *this;
483}
484
fe4da5cc 485
486