ATLAS Offline Software
Loading...
Searching...
No Matches
Pythia8::ZprimeFlatpT Class Reference
Inheritance diagram for Pythia8::ZprimeFlatpT:
Collaboration diagram for Pythia8::ZprimeFlatpT:

Public Member Functions

 ZprimeFlatpT ()
 ~ZprimeFlatpT ()
virtual bool canModifySigma ()
virtual double multiplySigmaBy (const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool)

Private Attributes

Pythia8_UserHooks::UserSetting< double > m_maxSHat = Pythia8_UserHooks::UserSetting<double>("ZprimeFlatpT:MaxSHat", 3000.)
Pythia8_UserHooks::UserSetting< double > m_doDecayWeightBelow = Pythia8_UserHooks::UserSetting<double>("ZprimeFlatpT:DoDecayWeightBelow", 0.)

Detailed Description

Definition at line 17 of file ZprimeFlatpT.cxx.

Constructor & Destructor Documentation

◆ ZprimeFlatpT()

Pythia8::ZprimeFlatpT::ZprimeFlatpT ( )
inline

Definition at line 22 of file ZprimeFlatpT.cxx.

22{}

◆ ~ZprimeFlatpT()

Pythia8::ZprimeFlatpT::~ZprimeFlatpT ( )
inline

Definition at line 25 of file ZprimeFlatpT.cxx.

25{}

Member Function Documentation

◆ canModifySigma()

virtual bool Pythia8::ZprimeFlatpT::canModifySigma ( )
inlinevirtual

Definition at line 28 of file ZprimeFlatpT.cxx.

28{return true;}

◆ multiplySigmaBy()

virtual double Pythia8::ZprimeFlatpT::multiplySigmaBy ( const SigmaProcess * sigmaProcessPtr,
const PhaseSpace * phaseSpacePtr,
bool  )
inlinevirtual

Definition at line 31 of file ZprimeFlatpT.cxx.

33 {
34 // All events should be 2 -> 1, but kill them if not.
35 // if (sigmaProcessPtr->nFinal() != 1) return 0.;
36 if (sigmaProcessPtr->nFinal() != 1) return 0.;
37
38 // Weight cross section with BW propagator, i.e. to remove it.
39 // (inEvent = false for initialization).
40 // No inEvent criteria, want weight both for cross section
41 // and MC generation.
42 int idRes = sigmaProcessPtr->resonanceB();
43 double mRes = particleDataPtr->m0(idRes);
44
45 double wRes = particleDataPtr->mWidth(idRes);
46 double m2Res = mRes*mRes;
47 double GamMRat = wRes/mRes;
48 double sHat = phaseSpacePtr->sHat();
49 long double weightBW = pow2(sHat - m2Res) + pow2(sHat * GamMRat);
50 long double weightpT= 1;
51 double rH = std::sqrt(sHat);
52
53 double par[2]={-8.95719e+00,1.62584e-03};
54 weightpT=std::exp(par[0]+par[1]*rH);
55 if(rH>=m_maxSHat(settingsPtr)){ weightpT=0; }
56
57 double weightDecay = 1.;
58
59 if(rH < m_doDecayWeightBelow(settingsPtr)){
60 double p0 = -0.000527117;
61 double p1 = 2.64665e-06;
62
63 weightDecay = 0.008/(p0+(p1*rH));
64 }
65
66 long double weight = weightBW * weightpT * weightDecay;
67 return weight;
68
69 }
static const std::map< unsigned int, unsigned int > pow2
Pythia8_UserHooks::UserSetting< double > m_maxSHat
Pythia8_UserHooks::UserSetting< double > m_doDecayWeightBelow

Member Data Documentation

◆ m_doDecayWeightBelow

Pythia8_UserHooks::UserSetting<double> Pythia8::ZprimeFlatpT::m_doDecayWeightBelow = Pythia8_UserHooks::UserSetting<double>("ZprimeFlatpT:DoDecayWeightBelow", 0.)
private

Definition at line 79 of file ZprimeFlatpT.cxx.

◆ m_maxSHat

Pythia8_UserHooks::UserSetting<double> Pythia8::ZprimeFlatpT::m_maxSHat = Pythia8_UserHooks::UserSetting<double>("ZprimeFlatpT:MaxSHat", 3000.)
private

Definition at line 75 of file ZprimeFlatpT.cxx.


The documentation for this class was generated from the following file: