CMS 3D CMS Logo

Py8PtotGun.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
7 
8 namespace gen {
9 
10  class Py8PtotGun : public Py8GunBase {
11  public:
13  ~Py8PtotGun() override {}
14 
15  bool generatePartonsAndHadronize() override;
16  const char* classname() const override;
17 
18  private:
19  // Ptot Gun particle(s) characteristics
20  double fMinEta;
21  double fMaxEta;
22  double fMinPtot;
23  double fMaxPtot;
25  };
26 
27  // implementation
28  //
30  // ParameterSet defpset ;
31  edm::ParameterSet pgun_params = ps.getParameter<edm::ParameterSet>("PGunParameters"); // , defpset ) ;
32  fMinEta = pgun_params.getParameter<double>("MinEta"); // ,-2.2);
33  fMaxEta = pgun_params.getParameter<double>("MaxEta"); // , 2.2);
34  fMinPtot = pgun_params.getParameter<double>("MinPtot"); // , 0.);
35  fMaxPtot = pgun_params.getParameter<double>("MaxPtot"); // , 0.);
36  fAddAntiParticle = pgun_params.getParameter<bool>("AddAntiParticle"); //, false) ;
37  }
38 
40  fMasterGen->event.reset();
41 
42  for (size_t i = 0; i < fPartIDs.size(); i++) {
43  int particleID = fPartIDs[i]; // this is PDG - need to convert to Py8 ???
44 
45  double phi = (fMaxPhi - fMinPhi) * randomEngine().flat() + fMinPhi;
46  double eta = (fMaxEta - fMinEta) * randomEngine().flat() + fMinEta;
47  double the = 2. * atan(exp(-eta));
48 
49  double pp = (fMaxPtot - fMinPtot) * randomEngine().flat() + fMinPtot;
50 
51  double mass = (fMasterGen->particleData).m0(particleID);
52 
53  double pt = pp * sin(the);
54  double ee = sqrt(pp * pp + mass * mass);
55  double px = pt * cos(phi);
56  double py = pt * sin(phi);
57  double pz = pp * cos(the);
58 
59  if (!((fMasterGen->particleData).isParticle(particleID))) {
61  }
62  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) // quarks
63  (fMasterGen->event).append(particleID, 23, 101, 0, px, py, pz, ee, mass);
64  else if (std::abs(particleID) == 21) // gluons
65  (fMasterGen->event).append(21, 23, 101, 102, px, py, pz, ee, mass);
66  // other
67  else {
68  (fMasterGen->event).append(particleID, 1, 0, 0, px, py, pz, ee, mass);
69  int eventSize = (fMasterGen->event).size() - 1;
70  // -log(flat) = exponential distribution
71  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
72  (fMasterGen->event)[eventSize].tau(tauTmp);
73  }
74 
75  // Here also need to add anti-particle (if any)
76  // otherwise just add a 2nd particle of the same type
77  // (for example, gamma)
78  //
79  if (fAddAntiParticle) {
80  if (1 <= std::abs(particleID) && std::abs(particleID) <= 6) { // quarks
81  (fMasterGen->event).append(-particleID, 23, 0, 101, -px, -py, -pz, ee, mass);
82  } else if (std::abs(particleID) == 21) { // gluons
83  (fMasterGen->event).append(21, 23, 102, 101, -px, -py, -pz, ee, mass);
84  } else {
85  if ((fMasterGen->particleData).isParticle(-particleID)) {
86  (fMasterGen->event).append(-particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
87  } else {
88  (fMasterGen->event).append(particleID, 1, 0, 0, -px, -py, -pz, ee, mass);
89  }
90  int eventSize = (fMasterGen->event).size() - 1;
91  // -log(flat) = exponential distribution
92  double tauTmp = -(fMasterGen->event)[eventSize].tau0() * log(randomEngine().flat());
93  (fMasterGen->event)[eventSize].tau(tauTmp);
94  }
95  }
96  }
97 
98  if (!fMasterGen->next())
99  return false;
100  evtGenDecay();
101 
102  event() = std::make_unique<HepMC::GenEvent>();
103  return toHepMC.fill_next_event(fMasterGen->event, event().get());
104  }
105 
106  const char* Py8PtotGun::classname() const { return "Py8PtotGun"; }
107 
109 
110 } // namespace gen
111 
112 using gen::Pythia8PtotGun;
gen::Py8InterfaceBase::fMasterGen
std::unique_ptr< Pythia8::Pythia > fMasterGen
Definition: Py8InterfaceBase.h:47
gen::Py8PtotGun::fMaxPtot
double fMaxPtot
Definition: Py8PtotGun.cc:23
mps_fire.i
i
Definition: mps_fire.py:428
gen::Py8GunBase::fMinPhi
double fMinPhi
Definition: Py8GunBase.h:58
gen::Py8PtotGun::fAddAntiParticle
bool fAddAntiParticle
Definition: Py8PtotGun.cc:24
metsig::tau
Definition: SignAlgoResolutions.h:49
gen::Py8PtotGun::classname
const char * classname() const override
Definition: Py8PtotGun.cc:106
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
gen::P8RndmEngine::flat
double flat() override
Definition: P8RndmEngine.cc:7
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
gen::Py8PtotGun::generatePartonsAndHadronize
bool generatePartonsAndHadronize() override
Definition: Py8PtotGun.cc:39
ExternalDecayDriver.h
gen::Py8PtotGun::fMinPtot
double fMinPtot
Definition: Py8PtotGun.cc:22
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
gen::Py8PtotGun
Definition: Py8PtotGun.cc:10
EgammaObjectsElectrons_cfi.particleID
particleID
Definition: EgammaObjectsElectrons_cfi.py:4
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
gen::Py8PtotGun::~Py8PtotGun
~Py8PtotGun() override
Definition: Py8PtotGun.cc:13
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::GeneratorFilter
Definition: GeneratorFilter.h:40
gen::Py8GunBase::evtGenDecay
void evtGenDecay()
Definition: Py8GunBase.cc:144
PVValHelper::eta
Definition: PVValidationHelpers.h:70
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
gen::Py8PtotGun::Py8PtotGun
Py8PtotGun(edm::ParameterSet const &)
Definition: Py8PtotGun.cc:29
gen
Definition: PythiaDecays.h:13
gen::Py8InterfaceBase::randomEngine
P8RndmEngine & randomEngine()
Definition: Py8InterfaceBase.h:44
gen::Py8PtotGun::fMaxEta
double fMaxEta
Definition: Py8PtotGun.cc:21
edm::ParameterSet
Definition: ParameterSet.h:47
GeneratorFilter.h
mps_setup.append
append
Definition: mps_setup.py:85
gen::Py8GunBase::fPartIDs
std::vector< int > fPartIDs
Definition: Py8GunBase.h:57
get
#define get
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
gen::Pythia8PtotGun
edm::GeneratorFilter< gen::Py8PtotGun, gen::ExternalDecayDriver > Pythia8PtotGun
Definition: Py8PtotGun.cc:108
gen::BaseHadronizer::event
std::unique_ptr< HepMC::GenEvent > & event()
Definition: BaseHadronizer.h:86
gen::Py8GunBase
Definition: Py8GunBase.h:40
gen::Py8InterfaceBase::toHepMC
HepMC::Pythia8ToHepMC toHepMC
Definition: Py8InterfaceBase.h:49
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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
gen::Py8PtotGun::fMinEta
double fMinEta
Definition: Py8PtotGun.cc:20
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443