CMS 3D CMS Logo

FlatRandomMultiParticlePGunProducer.cc
Go to the documentation of this file.
1 #include <ostream>
2 
4 
7 
13 #include "CLHEP/Random/RandFlat.h"
14 
15 //#define DebugLog
16 using namespace edm;
17 
19 
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()) throw cms::Exception("Configuration") << "Not all probabilities given for all particle types " << fProbParticle_.size() << ":" << fPartIDs.size() << " need them to match\n";
25 
26  produces<HepMCProduct>("unsmeared");
27  produces<GenEventInfoProduct>();
28 #ifdef DebugLog
29  std::cout << "Internal FlatRandomPGun is initialzed for " << fPartIDs.size()
30  << " particles in momentum range " << fMinP_ << ":" << fMaxP_
31  << std::endl ;
32  for (unsigned int k=0; k<fPartIDs.size(); ++k)
33  std::cout << " [" << k << "] " << fPartIDs[k] << ":" << fProbParticle_[k];
34  std::cout << std::endl;
35 #endif
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 #ifdef DebugLog
41  std::cout << "Corrected probabilities:";
42  for (unsigned int k=0; k<fProbParticle_.size(); ++k)
43  std::cout << " " << fProbParticle_[k];
44  std::cout << std::endl;
45 #endif
46 }
47 
49 
51  const edm::EventSetup& es) {
52 
54  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
55 
56 #ifdef DebugLog
57  if (fVerbosity > 0)
58  std::cout << "FlatRandomMultiParticlePGunProducer: Begin New Event Generation" << std::endl ;
59 #endif
60 
61  // event loop (well, another step in it...)
62  // no need to clean up GenEvent memory - done in HepMCProduct
63  // here re-create fEvt (memory)
64  //
65  fEvt = new HepMC::GenEvent() ;
66 
67  // now actualy, cook up the event from PDGTable and gun parameters
68  //
69 
70  // 1st, primary vertex
71  //
72  HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
73 
74  // loop over particles
75  //
76  int barcode(0), PartID(fPartIDs[0]);
77  double r1 = CLHEP::RandFlat::shoot(engine, 0., 1.);
78  for (unsigned int ip=0; ip<fPartIDs.size(); ip++) {
79  if (r1 <= fProbParticle_[ip]) {
80  PartID = fPartIDs[ip];
81  break;
82  }
83  }
84 #ifdef DebugLog
85  if (fVerbosity > 0) std::cout << "Random " << r1 << " PartID " << PartID
86  << std::endl;
87 #endif
88  double mom = CLHEP::RandFlat::shoot(engine, fMinP_, fMaxP_);
89  double eta = CLHEP::RandFlat::shoot(engine, fMinEta, fMaxEta);
90  double phi = CLHEP::RandFlat::shoot(engine, fMinPhi, fMaxPhi);
91  const HepPDT::ParticleData*
92  PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
93  double mass = PData->mass().value() ;
94  double energy = sqrt(mom*mom + mass*mass);
95  double theta = 2.*atan(exp(-eta));
96  double px = mom*sin(theta)*cos(phi);
97  double py = mom*sin(theta)*sin(phi);
98  double pz = mom*cos(theta);
99 
100  HepMC::FourVector p(px,py,pz,energy) ;
101  HepMC::GenParticle* Part = new HepMC::GenParticle(p,PartID,1);
102  barcode++ ;
103  Part->suggest_barcode(barcode);
104  Vtx->add_particle_out(Part);
105 
106  if (fAddAntiParticle) {
107  HepMC::FourVector ap(-px,-py,-pz,energy) ;
108  int APartID = (PartID == 22 || PartID == 23) ? PartID : -PartID;
109  HepMC::GenParticle* APart = new HepMC::GenParticle(ap,APartID,1);
110  barcode++ ;
111  APart->suggest_barcode(barcode);
112  Vtx->add_particle_out(APart);
113  }
114 
115  fEvt->add_vertex(Vtx) ;
116  fEvt->set_event_number(e.id().event()) ;
117  fEvt->set_signal_process_id(20) ;
118 
119 #ifdef DebugLog
120  if (fVerbosity > 0) fEvt->print() ;
121 #endif
122 
123  std::unique_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
124  BProduct->addHepMCData( fEvt );
125  e.put(std::move(BProduct), "unsmeared");
126 
127  std::unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
128  e.put(std::move(genEventInfo));
129 #ifdef DebugLog
130  if (fVerbosity > 0)
131  std::cout << " FlatRandomMultiParticlePGunProducer : Event Generation Done " << std::endl;
132 #endif
133 }
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
ESHandle< HepPDT::ParticleDataTable > fPDGTable
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > fPartIDs
HepPDT::ParticleData ParticleData
int k[5][pyjets_maxn]
edm::EventID id() const
Definition: EventBase.h:60
HLT enums.
StreamID streamID() const
Definition: Event.h:81
def move(src, dest)
Definition: eostools.py:510
virtual void produce(Event &e, const EventSetup &es) override