CMS 3D CMS Logo

RandomMultiParticlePGunProducer.cc
Go to the documentation of this file.
1 #include <ostream>
2 
4 
7 
14 #include "CLHEP/Random/RandFlat.h"
15 
16 //#define DebugLog
17 using namespace edm;
18 
20 
21  ParameterSet pgunParams = pset.getParameter<ParameterSet>("PGunParameters");
22  fProbParticle_ = pgunParams.getParameter<std::vector<double> >("ProbParts");
23  fProbP_ = pgunParams.getParameter<std::vector<double> >("ProbP");
24  fMinP_ = pgunParams.getParameter<double>("MinP");
25  fMaxP_ = pgunParams.getParameter<double>("MaxP");
26  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";
27 
28  produces<HepMCProduct>("unsmeared");
29  produces<GenEventInfoProduct>();
30  fBinsP_ = (int)(fProbP_.size());
31  fDeltaP_ = (fMaxP_-fMinP_)/fBinsP_;
32 #ifdef DebugLog
33  edm::LogVerbatim("IOMC") << "Internal FlatRandomPGun is initialzed for "
34  << fPartIDs.size() << " particles in "
35  << fBinsP_ << " bins within momentum range "
36  << fMinP_ << ":" << fMaxP_;
37  for (unsigned int k=0; k<fPartIDs.size(); ++k)
38  edm::LogVerbatim("IOMC") << " [" << k << "] " << fPartIDs[k] << ":"
39  << fProbParticle_[k];
40  edm::LogVerbatim("IOMC") << "Momentum distribution is given by";
41  for (int k=0; k<fBinsP_; ++k) {
42  double p = fMinP_ + k*fDeltaP_;
43  edm::LogVerbatim("IOMC") << " Bin[" << k << "] " << p << ":" << p+fDeltaP_
44  << " --> " << fProbP_[k];
45  }
46 #endif
47  for (unsigned int k=1; k<fProbParticle_.size(); ++k)
48  fProbParticle_[k] += fProbParticle_[k-1];
49  for (unsigned int k=0; k<fProbParticle_.size(); ++k)
50  fProbParticle_[k] /= fProbParticle_[fProbParticle_.size()-1];
51 #ifdef DebugLog
52  edm::LogVerbatim("IOMC") << "Corrected probabilities for particle type:";
53  for (unsigned int k=0; k<fProbParticle_.size(); ++k)
54  edm::LogVerbatim("IOMC") << " [" << k << "]: " << fProbParticle_[k];
55 #endif
56  for (int k=1; k<fBinsP_; ++k)
57  fProbP_[k] += fProbP_[k-1];
58  for (int k=0; k<fBinsP_; ++k)
59  fProbP_[k] /= fProbP_[fBinsP_-1];
60 #ifdef DebugLog
61  edm::LogVerbatim("IOMC") << "Corrected probabilities for momentum:";
62  for (int k=0; k<fBinsP_; ++k) {
63  double p = fMinP_ + k*fDeltaP_;
64  edm::LogVerbatim("IOMC") << " Bin[" << k << "] " << p << ":" << p+fDeltaP_
65  << " --> " << fProbP_[k];
66  }
67 #endif
68 }
69 
71  const edm::EventSetup& es) {
72 
74  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
75 
76 #ifdef DebugLog
77  if (fVerbosity > 0)
78  edm::LogVerbatim("IOMC") << "RandomMultiParticlePGunProducer: "
79  << "Begin New Event Generation";
80 #endif
81 
82  // event loop (well, another step in it...)
83  // no need to clean up GenEvent memory - done in HepMCProduct
84  // here re-create fEvt (memory)
85  //
86  fEvt = new HepMC::GenEvent() ;
87 
88  // now actualy, cook up the event from PDGTable and gun parameters
89  //
90 
91  // 1st, primary vertex
92  //
93  HepMC::GenVertex* Vtx = new HepMC::GenVertex( HepMC::FourVector(0.,0.,0.));
94 
95  // loop over particles
96  //
97  int barcode(0), PartID(fPartIDs[0]);
98  double r1 = CLHEP::RandFlat::shoot(engine, 0., 1.);
99  for (unsigned int ip=0; ip<fPartIDs.size(); ip++) {
100  if (r1 <= fProbParticle_[ip]) {
101  PartID = fPartIDs[ip];
102  break;
103  }
104  }
105  double minP(fMinP_);
106  double r2 = CLHEP::RandFlat::shoot(engine, 0., 1.);
107  for (int ip=0; ip<fBinsP_; ip++) {
108  if (r2 <= fProbP_[ip]) {
109  minP = fMinP_ + ip*fDeltaP_;
110  break;
111  }
112  }
113  double maxP = minP+fDeltaP_;
114 #ifdef DebugLog
115  if (fVerbosity > 0)
116  edm::LogVerbatim("IOMC") << "Random " << r1 << " PartID " << PartID
117  << " and for p " << r2 << " range " << minP
118  << ":" << maxP;
119 #endif
120  double mom = CLHEP::RandFlat::shoot(engine, minP, maxP);
121  double eta = CLHEP::RandFlat::shoot(engine, fMinEta, fMaxEta);
122  double phi = CLHEP::RandFlat::shoot(engine, fMinPhi, fMaxPhi);
123  const HepPDT::ParticleData*
124  PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
125  double mass = PData->mass().value() ;
126  double energy = sqrt(mom*mom + mass*mass);
127  double theta = 2.*atan(exp(-eta));
128  double px = mom*sin(theta)*cos(phi);
129  double py = mom*sin(theta)*sin(phi);
130  double pz = mom*cos(theta);
131 
132  HepMC::FourVector p(px,py,pz,energy) ;
133  HepMC::GenParticle* Part = new HepMC::GenParticle(p,PartID,1);
134  barcode++ ;
135  Part->suggest_barcode(barcode);
136  Vtx->add_particle_out(Part);
137 
138  if (fAddAntiParticle) {
139  HepMC::FourVector ap(-px,-py,-pz,energy) ;
140  int APartID = (PartID == 22 || PartID == 23) ? PartID : -PartID;
141  HepMC::GenParticle* APart = new HepMC::GenParticle(ap,APartID,1);
142  barcode++ ;
143  APart->suggest_barcode(barcode);
144  Vtx->add_particle_out(APart);
145  }
146 
147  fEvt->add_vertex(Vtx) ;
148  fEvt->set_event_number(e.id().event()) ;
149  fEvt->set_signal_process_id(20) ;
150 
151 #ifdef DebugLog
152  if (fVerbosity > 0) fEvt->print() ;
153 #endif
154 
155  std::unique_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
156  BProduct->addHepMCData( fEvt );
157  e.put(std::move(BProduct), "unsmeared");
158 
159  std::unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
160  e.put(std::move(genEventInfo));
161 #ifdef DebugLog
162  if (fVerbosity > 0)
163  edm::LogVerbatim("IOMC") << " RandomMultiParticlePGunProducer : Event Generation Done ";
164 #endif
165 }
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:137
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.
void produce(Event &e, const EventSetup &es) override
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:96
def move(src, dest)
Definition: eostools.py:510