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