BmnRoot
Loading...
Searching...
No Matches
MpdHypYPtGenerator.cxx
Go to the documentation of this file.
1// -------------------------------------------------------------------------
2// ----- CbmAnaHypYPtGenerator source file -----
3// ----- Created 03/10/04 by E. Kryshen -----
4// ----- Updated 10/02/10 by A. Zinchenko for MPD (with name change) -----
5// -------------------------------------------------------------------------
6
7#include "TRandom.h"
9#include "FairPrimaryGenerator.h"
10#include "TParticlePDG.h"
11#include "TDatabasePDG.h"
12#include "TF1.h"
13#include "TMath.h"
14
15// ------------------------------------------------------------------------
17{
18 // Default constructor
19
20 fPDGType = -1;
21 fMult = 0;
22}
23
24// ------------------------------------------------------------------------
25MpdHypYPtGenerator::MpdHypYPtGenerator(Int_t pdgid, Int_t mult, Double_t yield)
26 : FairGenerator(), fPDGType(pdgid), fMult(mult), fYield(yield)
27{
28 // Constructor. Set default distributions
29
32 SetRangePt();
33}
34// ------------------------------------------------------------------------
36{
37 // Initialize generator
38
39 // Check for particle type
40 TDatabasePDG *pdgBase = TDatabasePDG::Instance();
41 TParticlePDG *particle = pdgBase->GetParticle(fPDGType);
42 if (!particle) {
43 Fatal("CbmAnaHypYPtGenerator", "PDG code %d not defined.", fPDGType);
44 return kFALSE;
45 }
46 fPDGMass = particle->Mass();
47 // gRandom->SetSeed(0);
48 fDistPt = new TF1("distPt", "x*exp(-sqrt(x*x+[1]*[1])/[0])", fPtMin, fPtMax);
49 fDistPt->SetParameters(fT, fPDGMass, fY0, fSigma);
50 Info("Init", "pdg=%i y0=%4.2f sigma_y=%4.2f T_pt=%6.4f", fPDGType, fY0, fSigma, fT);
51 return kTRUE;
52}
53
54// ------------------------------------------------------------------------
55Bool_t MpdHypYPtGenerator::ReadEvent(FairPrimaryGenerator *primGen)
56{
57
58 Double_t phi, pt, y, mt, px, py, pz;
59 gRandom->Rndm();
60 if (fYield > 0 && gRandom->Rndm() > fYield) return kTRUE;
61
62 // Generate particles
63
64 for (Int_t k = 0; k < fMult; k++) {
65
66 phi = gRandom->Uniform(0, TMath::TwoPi());
67 pt = fDistPt->GetRandom(fPtMin, fPtMax);
68 px = pt * TMath::Cos(phi);
69 py = pt * TMath::Sin(phi);
70 y = gRandom->Gaus(fY0, fSigma);
71 mt = TMath::Sqrt(fPDGMass * fPDGMass + pt * pt);
72 pz = mt * TMath::SinH(y);
73 Info("ReadEvent", "Particle generated: pdg=%i pt=%f y=%f\n", fPDGType, pt, y);
74 primGen->AddTrack(fPDGType, px, py, pz, 0, 0, 0);
75 }
76 return kTRUE;
77}
78// ------------------------------------------------------------------------
virtual Bool_t ReadEvent(FairPrimaryGenerator *primGen)
void SetDistributionPt(Double_t T=0.223)
void SetRangePt(Double_t ptMin=0, Double_t ptMax=3)
void SetDistributionY(Double_t y0=0, Double_t sigma=0.72)