CMS 3D CMS Logo

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