CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FlatRandomMultiParticlePGunProducer.cc
Go to the documentation of this file.
1 #include <ostream>
2 
4 
7 
14 #include "CLHEP/Random/RandFlat.h"
15 
16 using namespace edm;
17 
19  : BaseFlatGunProducer(pset) {
20  ParameterSet pgunParams = pset.getParameter<ParameterSet>("PGunParameters");
21  fProbParticle_ = pgunParams.getParameter<std::vector<double> >("ProbParts");
22  fMinP_ = pgunParams.getParameter<double>("MinP");
23  fMaxP_ = pgunParams.getParameter<double>("MaxP");
24  if (fProbParticle_.size() != fPartIDs.size())
25  throw cms::Exception("Configuration") << "Not all probabilities given for all particle types "
26  << fProbParticle_.size() << ":" << fPartIDs.size() << " need them to match\n";
27 
28  produces<HepMCProduct>("unsmeared");
29  produces<GenEventInfoProduct>();
30 
31  edm::LogVerbatim("ParticleGun") << "FlatRandomMultiParticlePGun is initialzed for " << fPartIDs.size()
32  << " particles in momentum range " << fMinP_ << ":" << fMaxP_;
33  for (unsigned int k = 0; k < fPartIDs.size(); ++k)
34  edm::LogVerbatim("ParticleGun") << " [" << k << "] " << fPartIDs[k] << ":" << fProbParticle_[k];
35 
36  for (unsigned int k = 1; k < fProbParticle_.size(); ++k)
37  fProbParticle_[k] += fProbParticle_[k - 1];
38  for (unsigned int k = 0; k < fProbParticle_.size(); ++k) {
39  fProbParticle_[k] /= fProbParticle_[fProbParticle_.size() - 1];
40  edm::LogVerbatim("ParticleGun") << "Corrected probaility for " << fPartIDs[k] << ":" << fProbParticle_[k];
41  }
42 }
43 
45 
48  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
49 
50  if (fVerbosity > 0)
51  edm::LogVerbatim("ParticleGun") << "FlatRandomMultiParticlePGunProducer: Begin New Event Generation";
52 
53  // event loop (well, another step in it...)
54  // no need to clean up GenEvent memory - done in HepMCProduct
55  // here re-create fEvt (memory)
56  //
57  fEvt = new HepMC::GenEvent();
58 
59  // now actualy, cook up the event from PDGTable and gun parameters
60  //
61 
62  // 1st, primary vertex
63  //
64  HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0., 0., 0.));
65 
66  // loop over particles
67  //
68  int barcode(0), PartID(fPartIDs[0]);
69  double r1 = CLHEP::RandFlat::shoot(engine, 0., 1.);
70  for (unsigned int ip = 0; ip < fPartIDs.size(); ip++) {
71  if (r1 <= fProbParticle_[ip])
72  break;
73  PartID = fPartIDs[ip];
74  }
75  if (fVerbosity > 0)
76  edm::LogVerbatim("ParticleGun") << "Random " << r1 << " PartID " << PartID;
77 
78  double mom = CLHEP::RandFlat::shoot(engine, fMinP_, fMaxP_);
79  double eta = CLHEP::RandFlat::shoot(engine, fMinEta, fMaxEta);
80  double phi = CLHEP::RandFlat::shoot(engine, fMinPhi, fMaxPhi);
81  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID)));
82  double mass = PData->mass().value();
83  double energy = sqrt(mom * mom + mass * mass);
84  double theta = 2. * atan(exp(-eta));
85  double px = mom * sin(theta) * cos(phi);
86  double py = mom * sin(theta) * sin(phi);
87  double pz = mom * cos(theta);
88 
89  HepMC::FourVector p(px, py, pz, energy);
90  HepMC::GenParticle* Part = new HepMC::GenParticle(p, PartID, 1);
91  barcode++;
92  Part->suggest_barcode(barcode);
93  Vtx->add_particle_out(Part);
94 
95  if (fAddAntiParticle) {
96  HepMC::FourVector ap(-px, -py, -pz, energy);
97  int APartID = (PartID == 22 || PartID == 23) ? PartID : -PartID;
98  HepMC::GenParticle* APart = new HepMC::GenParticle(ap, APartID, 1);
99  barcode++;
100  APart->suggest_barcode(barcode);
101  Vtx->add_particle_out(APart);
102  }
103 
104  fEvt->add_vertex(Vtx);
105  fEvt->set_event_number(e.id().event());
106  fEvt->set_signal_process_id(20);
107 
108  if (fVerbosity > 1)
109  fEvt->print();
110 
111  std::unique_ptr<HepMCProduct> BProduct(new HepMCProduct());
112  BProduct->addHepMCData(fEvt);
113  e.put(std::move(BProduct), "unsmeared");
114 
115  std::unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
116  e.put(std::move(genEventInfo));
117 
118  if (fVerbosity > 0)
119  edm::LogVerbatim("ParticleGun") << " FlatRandomMultiParticlePGunProducer : Event Generation Done";
120 }
Log< level::Info, true > LogVerbatim
EventNumber_t event() const
Definition: EventID.h:40
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
ESHandle< HepPDT::ParticleDataTable > fPDGTable
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > fPartIDs
HepPDT::ParticleData ParticleData
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EventID id() const
Definition: EventBase.h:59
StreamID streamID() const
Definition: Event.h:98
void produce(Event &e, const EventSetup &es) override