1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
25 #include <TClonesArray.h>
26 #include <TDatabasePDG.h>
29 #include <TLorentzVector.h>
31 #include <TParticle.h>
32 #include <TParticlePDG.h>
34 #include <TVirtualMC.h>
36 #include "AliDecayer.h"
37 #include "AliGenMUONlib.h"
38 #include "AliGenParam.h"
41 #include "AliGenEventHeader.h"
45 //------------------------------------------------------------
49 <img src="picts/AliGenParam.gif">
53 //____________________________________________________________
54 AliGenParam::AliGenParam()
74 // Default constructor
76 //____________________________________________________________
77 AliGenParam::AliGenParam(Int_t npart, const AliGenLib * Library, Int_t param, const char* tname)
79 fPtParaFunc(Library->GetPt(param, tname)),
80 fYParaFunc (Library->GetY (param, tname)),
81 fIpParaFunc(Library->GetIp(param, tname)),
82 fV2ParaFunc(Library->GetV2(param, tname)),
98 // Constructor using number of particles parameterisation id and library
100 fTitle= "Particle Generator using pT and y parameterisation";
104 //____________________________________________________________
105 AliGenParam::AliGenParam(Int_t npart, Int_t param, const char* tname, const char* name):
126 // Constructor using parameterisation id and number of particles
129 fTitle= "Particle Generator using pT and y parameterisation";
131 AliGenLib* pLibrary = new AliGenMUONlib();
132 fPtParaFunc = pLibrary->GetPt(param, tname);
133 fYParaFunc = pLibrary->GetY (param, tname);
134 fIpParaFunc = pLibrary->GetIp(param, tname);
135 fV2ParaFunc = pLibrary->GetV2(param, tname);
139 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
142 SetChildMomentumRange();
145 SetChildThetaRange();
147 //____________________________________________________________
149 AliGenParam::AliGenParam(Int_t npart, Int_t param,
150 Double_t (*PtPara) (const Double_t*, const Double_t*),
151 Double_t (*YPara ) (const Double_t* ,const Double_t*),
152 Double_t (*V2Para) (const Double_t* ,const Double_t*),
153 Int_t (*IpPara) (TRandom *))
176 // Gines Martinez 1/10/99
178 fTitle= "Particle Generator using pT and y parameterisation";
182 for (Int_t i=0; i<5; i++) fChildSelect[i]=0;
185 SetChildMomentumRange();
188 SetChildThetaRange();
191 //____________________________________________________________
192 AliGenParam::~AliGenParam()
201 //-------------------------------------------------------------------
202 TVector3 AliGenParam::OrthogonalVector(TVector3 &inVec){
203 double abc[]={inVec.x(), inVec.y(), inVec.z()};
204 double xyz[]={1,1,1};
207 for(int i=0; i<3; i++)
212 xyz[solvDim]=(-abc[(1+solvDim)%3]-abc[(2+solvDim)%3])/abc[(0+solvDim)%3];
214 TVector3 res(xyz[0],xyz[1],xyz[2]);
218 void AliGenParam::RotateVector(Double_t *pin, Double_t *pout, Double_t costheta, Double_t sintheta,
219 Double_t cosphi, Double_t sinphi)
222 pout[0] = pin[0]*costheta*cosphi-pin[1]*sinphi+pin[2]*sintheta*cosphi;
223 pout[1] = pin[0]*costheta*sinphi+pin[1]*cosphi+pin[2]*sintheta*sinphi;
224 pout[2] = -1.0 * pin[0] * sintheta + pin[2] * costheta;
228 Double_t AliGenParam::IntegratedKrollWada(Double_t mh){
229 if(mh<0.003) return 0;
230 return 2*log(mh/0.000511/exp(7.0/4.0))/411.0/TMath::Pi();
233 double AliGenParam::ScreenFunction1(double screenVariable){
235 return 42.24 - 8.368 * log(screenVariable + 0.952);
237 return 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
240 double AliGenParam::ScreenFunction2(double screenVariable){
242 return 42.24 - 8.368 * log(screenVariable + 0.952);
244 return 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
247 double AliGenParam::RandomEnergyFraction(double Z, double photonEnergy){
250 double epsilon0Local = 0.000511 / photonEnergy ;
252 // Do it fast if photon energy < 2. MeV
253 if (photonEnergy < 0.002 )
255 epsilon = epsilon0Local + (0.5 - epsilon0Local) * fRandom->Rndm();
259 double fZ = 8*log(Z)/3;
260 double fcZ=(aZ*aZ)*(1/(1+aZ*aZ)+0.20206-0.0368*aZ*aZ+0.0083*aZ*aZ*aZ);
261 if (photonEnergy > 0.050) fZ += 8*fcZ;
263 // Limits of the screening variable
264 double screenFactor = 136. * epsilon0Local / pow (Z,1/3);
265 double screenMax = exp ((42.24 - fZ)/8.368) - 0.952 ;
266 double screenMin = std::min(4.*screenFactor,screenMax) ;
268 // Limits of the energy sampling
269 double epsilon1 = 0.5 - 0.5 * sqrt(1. - screenMin / screenMax) ;
270 double epsilonMin = std::max(epsilon0Local,epsilon1);
271 double epsilonRange = 0.5 - epsilonMin ;
273 // Sample the energy rate of the created electron (or positron)
277 double f10 = ScreenFunction1(screenMin) - fZ;
278 double f20 = ScreenFunction2(screenMin) - fZ;
279 double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
280 double normF2 = std::max(1.5 * f20,0.);
284 if (normF1 / (normF1 + normF2) > fRandom->Rndm() )
286 epsilon = 0.5 - epsilonRange * pow(fRandom->Rndm(), 0.333333) ;
287 screen = screenFactor / (epsilon * (1. - epsilon));
288 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
292 epsilon = epsilonMin + epsilonRange * fRandom->Rndm();
293 screen = screenFactor / (epsilon * (1 - epsilon));
294 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
296 } while ( gReject < fRandom->Rndm() );
298 } // End of epsilon sampling
302 double AliGenParam::RandomPolarAngle(){
304 const double a1 = 0.625;
308 // if (9. / (9. + d) > fRandom->Rndm())
309 if (0.25 > fRandom->Rndm())
311 u = - log(fRandom->Rndm() * fRandom->Rndm()) / a1 ;
315 u = - log(fRandom->Rndm() * fRandom->Rndm()) / a2 ;
320 Double_t AliGenParam::RandomMass(Double_t mh){
322 double y=fRandom->Rndm();
323 double mee=2*0.000511*TMath::Power(2*0.000511/mh,-y); //inverse of the enveloping cumulative distribution
324 double apxkw=2.0/3.0/137.036/TMath::Pi()/mee;
325 double val=fRandom->Uniform(0,apxkw); //enveloping probability density0
326 double kw=apxkw*sqrt(1-4*0.000511*0.000511/mee/mee)*(1+2*0.000511*0.000511/mee/mee)*1*1*TMath::Power(1-mee*mee/mh/mh,3);
332 Int_t AliGenParam::VirtualGammaPairProduction(TClonesArray *particles, Int_t nPart)
334 Int_t nPartNew=nPart;
335 for(int iPart=0; iPart<nPart; iPart++){
336 TParticle *gamma = (TParticle *) particles->At(iPart);
337 if(gamma->GetPdgCode()!=223000 || gamma->GetPdgCode()!=224000) continue;
338 //if(gamma->Energy()<0.001022) continue; //can never be
339 double mass=RandomMass(gamma->Pt());
341 // lepton pair kinematics in virtual photon rest frame
343 double Pe=TMath::Sqrt((Ee+0.000511)*(Ee-0.000511));
345 double costheta = (2.0 * gRandom->Rndm()) - 1.;
346 double sintheta = TMath::Sqrt((1. + costheta) * (1. - costheta));
347 double phi = 2.0 * TMath::ACos(-1.) * gRandom->Rndm();
348 double sinphi = TMath::Sin(phi);
349 double cosphi = TMath::Cos(phi);
351 // momentum vectors of leptons in virtual photon rest frame
352 Double_t pProd1[3] = {Pe * sintheta * cosphi,
353 Pe * sintheta * sinphi,
356 Double_t pProd2[3] = {-1.0 * Pe * sintheta * cosphi,
357 -1.0 * Pe * sintheta * sinphi,
358 -1.0 * Pe * costheta};
360 // lepton 4-vectors in properly rotated virtual photon rest frame
361 Double_t pRot1[3] = {0.};
362 RotateVector(pProd1, pRot1, costheta, -sintheta, -cosphi, -sinphi);
363 Double_t pRot2[3] = {0.};
364 RotateVector(pProd2, pRot2, costheta, -sintheta, -cosphi, -sinphi);
366 TLorentzVector e1V4(pRot1[0],pRot1[1],pRot1[2],Ee);
367 TLorentzVector e2V4(pRot2[0],pRot2[1],pRot2[2],Ee);
369 TVector3 boost(gamma->Px(),gamma->Py(),gamma->Pz());
370 boost*=1/sqrt(gamma->P()*gamma->P()+mass*mass);
375 gamma->ProductionVertex(vtx);
376 new((*particles)[nPartNew]) TParticle(11, gamma->GetStatusCode(), iPart+1, -1, 0, 0, e1V4, vtx);
378 new((*particles)[nPartNew]) TParticle(-11, gamma->GetStatusCode(), iPart+1, -1, 0, 0, e2V4, vtx);
383 Int_t AliGenParam::ForceGammaConversion(TClonesArray *particles, Int_t nPart)
385 //based on: http://geant4.cern.ch/G4UsersDocuments/UsersGuides/PhysicsReferenceManual/html/node27.html
386 // and: http://geant4.cern.ch/G4UsersDocuments/UsersGuides/PhysicsReferenceManual/html/node58.html
387 // and: G4LivermoreGammaConversionModel.cc
388 Int_t nPartNew=nPart;
389 for(int iPart=0; iPart<nPart; iPart++){
390 TParticle *gamma = (TParticle *) particles->At(iPart);
391 if(gamma->GetPdgCode()!=22) continue;
392 if(gamma->Energy()<0.001022) continue;
394 TVector3 gammaV3(gamma->Px(),gamma->Py(),gamma->Pz());
395 double frac=RandomEnergyFraction(1,gamma->Energy());
396 double Ee1=frac*gamma->Energy();
397 double Ee2=(1-frac)*gamma->Energy();
398 double Pe1=sqrt((Ee1+0.000511)*(Ee1-0.000511));
399 double Pe2=sqrt((Ee2+0.000511)*(Ee2-0.000511));
401 TVector3 rotAxis(OrthogonalVector(gammaV3));
402 Float_t az=fRandom->Uniform(TMath::Pi()*2);
403 rotAxis.Rotate(az,gammaV3);
404 TVector3 e1V3(gammaV3);
405 double u=RandomPolarAngle();
406 e1V3.Rotate(u/Ee1,rotAxis);
409 TVector3 e2V3(gammaV3);
410 e2V3.Rotate(-u/Ee2,rotAxis);
413 // gamma = new TParticle(*gamma);
414 // particles->RemoveAt(iPart);
415 gamma->SetFirstDaughter(nPartNew+1);
416 gamma->SetLastDaughter(nPartNew+2);
417 // new((*particles)[iPart]) TParticle(*gamma);
421 gamma->ProductionVertex(vtx);
422 new((*particles)[nPartNew]) TParticle(11, gamma->GetStatusCode(), iPart+1, -1, 0, 0, TLorentzVector(e1V3,Ee1), vtx);
424 new((*particles)[nPartNew]) TParticle(-11, gamma->GetStatusCode(), iPart+1, -1, 0, 0, TLorentzVector(e2V3,Ee2), vtx);
427 // particles->Compress();
428 return particles->GetEntriesFast();
431 //____________________________________________________________
432 void AliGenParam::Init()
436 if (gMC) fDecayer = gMC->GetDecayer();
439 <img src="picts/AliGenParam.gif">
443 snprintf(name, 256, "pt-parameterisation for %s", GetName());
445 if (fPtPara) fPtPara->Delete();
446 fPtPara = new TF1(name, fPtParaFunc, fPtMin, fPtMax,0);
447 gROOT->GetListOfFunctions()->Remove(fPtPara);
448 // Set representation precision to 10 MeV
449 Int_t npx= Int_t((fPtMax - fPtMin) / fDeltaPt);
451 fPtPara->SetNpx(npx);
453 snprintf(name, 256, "y-parameterisation for %s", GetName());
454 if (fYPara) fYPara->Delete();
455 fYPara = new TF1(name, fYParaFunc, fYMin, fYMax, 0);
456 gROOT->GetListOfFunctions()->Remove(fYPara);
458 snprintf(name, 256, "v2-parameterisation for %s", GetName());
459 if (fV2Para) fV2Para->Delete();
460 fV2Para = new TF1(name, fV2ParaFunc, fPtMin, fPtMax, 0);
461 // fV2Para = new TF1(name, "2*[0]/(1+TMath::Exp([1]*([2]-x)))-[0]", fPtMin, fPtMax);
462 // fV2Para->SetParameter(0, 0.236910);
463 // fV2Para->SetParameter(1, 1.71122);
464 // fV2Para->SetParameter(2, 0.0827617);
465 //gROOT->GetListOfFunctions()->Remove(fV2Para); //TR: necessary?
467 snprintf(name, 256, "dNdPhi for %s", GetName());
468 if (fdNdPhi) fdNdPhi->Delete();
469 fdNdPhi = new TF1(name, "1+2*[0]*TMath::Cos(2*(x-[1]))", fPhiMin, fPhiMax);
470 //gROOT->GetListOfFunctions()->Remove(fdNdPhi); //TR: necessary?
472 snprintf(name, 256, "pt-for-%s", GetName());
473 TF1 ptPara(name ,fPtParaFunc, 0, 15, 0);
474 snprintf(name, 256, "y-for-%s", GetName());
475 TF1 yPara(name, fYParaFunc, -6, 6, 0);
482 fdNdy0=fYParaFunc(&y1,&y2);
484 // Integral over generation region
485 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
486 Float_t intYS = yPara.Integral(fYMin, fYMax,(Double_t*) 0x0,1.e-6);
487 Float_t intPt0 = ptPara.Integral(0,15,(Double_t *) 0x0,1.e-6);
488 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,(Double_t*) 0x0,1.e-6);
490 Float_t intYS = yPara.Integral(fYMin, fYMax,1.e-6);
491 Float_t intPt0 = ptPara.Integral(0,15,1.e-6);
492 Float_t intPtS = ptPara.Integral(fPtMin,fPtMax,1.e-6);
494 Float_t phiWgt=(fPhiMax-fPhiMin)/2./TMath::Pi(); //TR: should probably be done differently in case of anisotropic phi...
495 if (fAnalog == kAnalog) {
496 fYWgt = intYS/fdNdy0;
497 fPtWgt = intPtS/intPt0;
498 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
500 fYWgt = intYS/fdNdy0;
501 fPtWgt = (fPtMax-fPtMin)/intPt0;
502 fParentWeight = fYWgt*fPtWgt*phiWgt/fNpart;
505 // particle decay related initialization
506 fDecayer->SetForceDecay(fForceDecay);
513 //____________________________________________________________
514 void AliGenParam::Generate()
517 // Generate 1 event (see Generate(Int_t ntimes) for details
521 //____________________________________________________________
522 void AliGenParam::GenerateN(Int_t ntimes)
525 // Generate ntimes*'npart' light and heavy mesons (J/Psi, upsilon or phi, Pion,
526 // Kaons, Etas, Omegas) and Baryons (proton, antiprotons, neutrons and
527 // antineutrons in the the desired theta, phi and momentum windows;
528 // Gaussian smearing on the vertex is done if selected.
529 // The decay of heavy mesons is done using lujet,
530 // and the childern particle are tracked by GEANT
531 // However, light mesons are directly tracked by GEANT
532 // setting fForceDecay = nodecay (SetForceDecay(nodecay))
535 // Reinitialize decayer
536 fDecayer->SetForceDecay(fForceDecay);
540 Float_t polar[3]= {0,0,0}; // Polarisation of the parent particle (for GEANT tracking)
541 Float_t origin0[3]; // Origin of the generated parent particle (for GEANT tracking)
542 Float_t time0; // Time0 of the generated parent particle
543 Float_t pt, pl, ptot; // Transverse, logitudinal and total momenta of the parent particle
544 Float_t phi, theta; // Phi and theta spherical angles of the parent particle momentum
546 och[3]; // Momentum, polarisation and origin of the children particles from lujet
551 static TClonesArray *particles;
553 if(!particles) particles = new TClonesArray("TParticle",1000);
555 TDatabasePDG *pDataBase = TDatabasePDG::Instance();
559 // Calculating vertex position per event
560 for (j=0;j<3;j++) origin0[j]=fOrigin[j];
562 if(fVertexSmear==kPerEvent) {
564 for (j=0;j<3;j++) origin0[j]=fVertex[j];
570 // Generating fNpart particles
573 Int_t nGen = fNpart*ntimes;
578 Int_t iPart = fIpParaFunc(fRandom);
581 // custom pdg codes for to destinguish direct photons
582 if((iPart>=221000) && (iPart<=229000)) iPart=22;
584 fChildWeight=(fDecayer->GetPartialBranchingRatio(iPart))*fParentWeight;
585 TParticlePDG *particle = pDataBase->GetParticle(iPart);
586 Float_t am = particle->Mass();
591 ty = TMath::TanH(fYPara->GetRandom());
595 if (fAnalog == kAnalog) {
596 pt=fPtPara->GetRandom();
600 pt=fPtMin+random[1]*(fPtMax-fPtMin);
602 wgtp=fParentWeight*fPtParaFunc(& ptd, &dummy);
603 wgtch=fChildWeight*fPtParaFunc(& ptd, &dummy);
605 xmt=sqrt(pt*pt+am*am);
606 if (TMath::Abs(ty)==1.) {
609 "Division by 0: Please check you rapidity range !");
614 //phi=fEvPlane; //align first particle of each event with event plane
616 double v2 = fV2Para->Eval(pt);
617 fdNdPhi->SetParameter(0,v2);
618 fdNdPhi->SetParameter(1,fEvPlane);
619 phi=fdNdPhi->GetRandom();
622 pl=xmt*ty/sqrt((1.-ty)*(1.+ty));
623 theta=TMath::ATan2(pt,pl);
625 if(theta<fThetaMin || theta>fThetaMax) continue;
626 ptot=TMath::Sqrt(pt*pt+pl*pl);
628 if(ptot<fPMin || ptot>fPMax) continue;
630 p[0]=pt*TMath::Cos(phi);
631 p[1]=pt*TMath::Sin(phi);
634 if(fVertexSmear==kPerTrack) {
638 fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
639 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
642 time0 = fTimeOrigin + fOsigma[2]/TMath::Ccgs()*
643 TMath::Cos(2*random[0]*TMath::Pi())*
644 TMath::Sqrt(-2*TMath::Log(random[1]));
647 // Looking at fForceDecay :
648 // if fForceDecay != none Primary particle decays using
649 // AliPythia and children are tracked by GEANT
651 // if fForceDecay == none Primary particle is tracked by GEANT
652 // (In the latest, make sure that GEANT actually does all the decays you want)
654 Bool_t decayed = kFALSE;
657 if (fForceDecay != kNoDecay) {
658 // Using lujet to decay particle
659 Float_t energy=TMath::Sqrt(ptot*ptot+am*am);
660 TLorentzVector pmom(p[0], p[1], p[2], energy);
661 fDecayer->Decay(iPart,&pmom);
663 // select decay particles
664 Int_t np=fDecayer->ImportParticles(particles);
667 if(fForceConv) np=ForceGammaConversion(particles,np);
668 if(iPart==223000 || iPart==224000){
669 // wgtp*=IntegratedKrollWada(pt);
670 // wgtch*=IntegratedKrollWada(pt);
671 // np=VirtualGammaPairProduction(particles,np)
675 // Selecting GeometryAcceptance for particles fPdgCodeParticleforAcceptanceCut;
676 if (fGeometryAcceptance)
677 if (!CheckAcceptanceGeometry(np,particles)) continue;
679 Int_t* pFlag = new Int_t[np];
680 Int_t* pParent = new Int_t[np];
681 Int_t* pSelected = new Int_t[np];
682 Int_t* trackIt = new Int_t[np];
684 for (i=0; i<np; i++) {
692 TParticle* iparticle = 0;
694 for (i = 1; i<np ; i++) {
696 iparticle = (TParticle *) particles->At(i);
697 Int_t kf = iparticle->GetPdgCode();
698 Int_t ks = iparticle->GetStatusCode();
702 ipF = iparticle->GetFirstDaughter();
703 ipL = iparticle->GetLastDaughter();
704 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
708 // flag decay products of particles with long life-time (c tau > .3 mum)
711 // TParticlePDG *particle = pDataBase->GetParticle(kf);
713 Double_t lifeTime = fDecayer->GetLifetime(kf);
714 // Double_t mass = particle->Mass();
715 // Double_t width = particle->Width();
716 if (lifeTime > (Double_t) fMaxLifeTime) {
717 ipF = iparticle->GetFirstDaughter();
718 ipL = iparticle->GetLastDaughter();
719 if (ipF > 0) for (j=ipF-1; j<ipL; j++) pFlag[j]=1;
728 if ((ChildSelected(TMath::Abs(kf)) || fForceDecay == kAll || fSelectAll) && trackIt[i])
731 pc[0]=iparticle->Px();
732 pc[1]=iparticle->Py();
733 pc[2]=iparticle->Pz();
734 Bool_t childok = KinematicSelection(iparticle, 1);
745 } // if child selection
747 } // decay particle loop
748 } // if decay products
751 if ((fCutOnChild && ncsel >0) || !fCutOnChild){
757 PushTrack(0, -1, iPart, p, origin0, polar, time0, kPPrimary, nt, wgtp, ((decayed)? 11 : 1));
765 for (i = 1; i < np; i++) {
767 TParticle* iparticle = (TParticle *) particles->At(i);
768 Int_t kf = iparticle->GetPdgCode();
769 Int_t ksc = iparticle->GetStatusCode();
770 Int_t jpa = iparticle->GetFirstMother()-1;
772 och[0] = origin0[0]+iparticle->Vx();
773 och[1] = origin0[1]+iparticle->Vy();
774 och[2] = origin0[2]+iparticle->Vz();
775 pc[0] = iparticle->Px();
776 pc[1] = iparticle->Py();
777 pc[2] = iparticle->Pz();
780 iparent = pParent[jpa];
785 PushTrack(fTrackIt * trackIt[i], iparent, kf,
787 time0 + iparticle->T(), kPDecay, nt, wgtch, ksc);
795 if (pFlag) delete[] pFlag;
796 if (pParent) delete[] pParent;
797 if (pSelected) delete[] pSelected;
798 if (trackIt) delete[] trackIt;
799 } // kinematic selection
800 else // nodecay option, so parent will be tracked by GEANT (pions, kaons, eta, omegas, baryons)
803 PushTrack(fTrackIt,-1,iPart,p,origin0,polar,time0,kPPrimary,nt,wgtp, 1);
811 SetHighWaterMark(nt);
813 AliGenEventHeader* header = new AliGenEventHeader("PARAM");
814 header->SetPrimaryVertex(fVertex);
815 header->SetInteractionTime(fTime);
816 header->SetNProduced(fNprimaries);
819 //____________________________________________________________________________________
820 Float_t AliGenParam::GetRelativeArea(Float_t ptMin, Float_t ptMax, Float_t yMin, Float_t yMax, Float_t phiMin, Float_t phiMax)
823 // Normalisation for selected kinematic region
825 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
827 fPtPara->Integral(ptMin,ptMax,(Double_t *)0,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),(Double_t *)0,1.e-6) *
828 fYPara->Integral(yMin,yMax,(Double_t *)0,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),(Double_t *)0,1.e-6) *
829 (phiMax-phiMin)/360.;
832 fPtPara->Integral(ptMin,ptMax,1.e-6) / fPtPara->Integral( fPtPara->GetXmin(), fPtPara->GetXmax(),1.e-6) *
833 fYPara->Integral(yMin,yMax,1.e-6)/fYPara->Integral(fYPara->GetXmin(),fYPara->GetXmax(),1.e-6) *
834 (phiMax-phiMin)/360.;
836 return TMath::Abs(ratio);
839 //____________________________________________________________________________________
841 void AliGenParam::Draw( const char * /*opt*/)
844 // Draw the pT and y Distributions
846 TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
850 fPtPara->GetHistogram()->SetXTitle("p_{T} (GeV)");
853 fYPara->GetHistogram()->SetXTitle("y");