CMS 3D CMS Logo

Py8MassGun.cc
Go to the documentation of this file.
1 
2 #include <memory>
3 
6 
8 
9 namespace gen {
10 
11  class Py8MassGun : public Py8GunBase {
12  public:
14  ~Py8MassGun() override {}
15 
16  bool generatePartonsAndHadronize() override;
17  const char* classname() const override;
18 
19  private:
20  // PtGun particle(s) characteristics
21  double fMinEta;
22  double fMaxEta;
23  double fMinP;
24  double fMaxP;
25  double fMinPt;
26  double fMaxPt;
27  double fMinM;
28  double fMaxM;
29  int fMomMode;
30  };
31 
32  // implementation
33  //
35  // ParameterSet defpset ;
36  edm::ParameterSet pgun_params = ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
37  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
38  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
39  fMinP = pgun_params.getParameter<double>("MinP"); // , 0.);
40  fMaxP = pgun_params.getParameter<double>("MaxP"); // , 0.);
41  fMinPt = pgun_params.getParameter<double>("MinPt"); // , 0.);
42  fMaxPt = pgun_params.getParameter<double>("MaxPt"); // , 0.);
43  fMinM = pgun_params.getParameter<double>("MinM"); // , 0.);
44  fMaxM = pgun_params.getParameter<double>("MaxM"); // , 0.);
45  fMomMode = pgun_params.getParameter<int>("MomMode"); // , 1);
46  }
47 
49  fMasterGen->event.reset();
50  size_t pSize = fPartIDs.size();
51  if (pSize > 2)
52  return false;
53 
54  // Pick a flat mass range
55  double phi, eta, the, ee, pp;
56  double m0 = (fMaxM - fMinM) * randomEngine().flat() + fMinM;
57  // Global eta
59 
60  if (pSize == 2) {
61  // Masses.
62  double m1 = fMasterGen->particleData.m0(fPartIDs[0]);
63  double m2 = fMasterGen->particleData.m0(fPartIDs[1]);
64 
65  // Energies and absolute momentum in the rest frame.
66  if (m1 + m2 > m0)
67  return false;
68  double e1 = 0.5 * (m0 * m0 + m1 * m1 - m2 * m2) / m0;
69  double e2 = 0.5 * (m0 * m0 + m2 * m2 - m1 * m1) / m0;
70  double pAbs = 0.5 * sqrt((m0 - m1 - m2) * (m0 + m1 + m2) * (m0 + m1 - m2) * (m0 - m1 + m2)) / m0;
71  // Isotropic angles in rest frame give three-momentum.
72  double cosTheta = 2. * randomEngine().flat() - 1.;
73  double sinTheta = sqrt(1. - cosTheta * cosTheta);
74  phi = 2. * M_PI * randomEngine().flat();
75 
76  double pX = pAbs * sinTheta * cos(phi);
77  double pY = pAbs * sinTheta * sin(phi);
78  double pZ = pAbs * cosTheta;
79 
80  (fMasterGen->event).append(fPartIDs[0], 1, 0, 0, pX, pY, pZ, e1, m1);
81  (fMasterGen->event).append(fPartIDs[1], 1, 0, 0, -pX, -pY, -pZ, e2, m2);
82  } else {
83  (fMasterGen->event).append(fPartIDs[0], 1, 0, 0, 0.0, 0.0, 0.0, m0, m0);
84  }
85 
86  //now the boost (from input params)
87  if (fMomMode == 0) {
88  pp = (fMaxP - fMinP) * randomEngine().flat() + fMinP;
89  } else {
90  double pT = (fMaxPt - fMinPt) * randomEngine().flat() + fMinPt;
91  pp = pT * cosh(eta);
92  }
93  ee = sqrt(m0 * m0 + pp * pp);
94 
95  //the boost direction (from input params)
96  //
97  phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
98  the = 2. * atan(exp(-eta));
99 
100  double betaX = pp / ee * std::sin(the) * std::cos(phi);
101  double betaY = pp / ee * std::sin(the) * std::sin(phi);
102  double betaZ = pp / ee * std::cos(the);
103 
104  // boost all particles
105  //
106  (fMasterGen->event).bst(betaX, betaY, betaZ);
107 
108  if (!fMasterGen->next())
109  return false;
110 
111  event() = std::make_unique<HepMC::GenEvent>();
112  return toHepMC.fill_next_event(fMasterGen->event, event().get());
113  }
114 
115  const char* Py8MassGun::classname() const { return "Py8MassGun"; }
116 
118 
119 } // namespace gen
120 
121 using gen::Pythia8MassGun;
gen::Py8InterfaceBase::fMasterGen
std::unique_ptr< Pythia8::Pythia > fMasterGen
Definition: Py8InterfaceBase.h:47
gen::Py8MassGun
Definition: Py8MassGun.cc:11
gen::Py8GunBase::fMinPhi
double fMinPhi
Definition: Py8GunBase.h:58
gen::Py8MassGun::classname
const char * classname() const override
Definition: Py8MassGun.cc:115
gen::P8RndmEngine::flat
double flat() override
Definition: P8RndmEngine.cc:7
gen::Py8MassGun::fMaxM
double fMaxM
Definition: Py8MassGun.cc:28
ExternalDecayDriver.h
gen::Py8MassGun::fMinPt
double fMinPt
Definition: Py8MassGun.cc:25
gen::Py8MassGun::~Py8MassGun
~Py8MassGun() override
Definition: Py8MassGun.cc:14
gen::Py8MassGun::Py8MassGun
Py8MassGun(edm::ParameterSet const &)
Definition: Py8MassGun.cc:34
gen::Py8MassGun::fMinP
double fMinP
Definition: Py8MassGun.cc:23
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
gen::Py8MassGun::generatePartonsAndHadronize
bool generatePartonsAndHadronize() override
Definition: Py8MassGun.cc:48
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::GeneratorFilter
Definition: GeneratorFilter.h:40
PVValHelper::eta
Definition: PVValidationHelpers.h:70
PVValHelper::pT
Definition: PVValidationHelpers.h:71
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
gen::Pythia8MassGun
edm::GeneratorFilter< gen::Py8MassGun, gen::ExternalDecayDriver > Pythia8MassGun
Definition: Py8MassGun.cc:117
gen
Definition: PythiaDecays.h:13
gen::Py8InterfaceBase::randomEngine
P8RndmEngine & randomEngine()
Definition: Py8InterfaceBase.h:44
edm::ParameterSet
Definition: ParameterSet.h:47
GeneratorFilter.h
mps_setup.append
append
Definition: mps_setup.py:85
StorageManager_cfg.e1
e1
Definition: StorageManager_cfg.py:16
gen::Py8MassGun::fMinEta
double fMinEta
Definition: Py8MassGun.cc:21
gen::Py8GunBase::fPartIDs
std::vector< int > fPartIDs
Definition: Py8GunBase.h:57
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
gen::Py8MassGun::fMaxEta
double fMaxEta
Definition: Py8MassGun.cc:22
get
#define get
gen::Py8MassGun::fMinM
double fMinM
Definition: Py8MassGun.cc:27
callgraph.m2
m2
Definition: callgraph.py:38
gen::BaseHadronizer::event
std::unique_ptr< HepMC::GenEvent > & event()
Definition: BaseHadronizer.h:86
gen::Py8MassGun::fMaxPt
double fMaxPt
Definition: Py8MassGun.cc:26
gen::Py8GunBase
Definition: Py8GunBase.h:40
gen::Py8InterfaceBase::toHepMC
HepMC::Pythia8ToHepMC toHepMC
Definition: Py8InterfaceBase.h:49
gen::Py8MassGun::fMomMode
int fMomMode
Definition: Py8MassGun.cc:29
gen::Py8MassGun::fMaxP
double fMaxP
Definition: Py8MassGun.cc:24
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
createTree.pp
pp
Definition: createTree.py:17
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
gen::Py8GunBase::fMaxPhi
double fMaxPhi
Definition: Py8GunBase.h:59
Py8GunBase.h
event
Definition: event.py:1