CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MultiParticleInConeGunProducer.cc
Go to the documentation of this file.
1 /*
2  * \author Jean-Roch Vlimant
3  */
4 
5 #include <ostream>
6 
8 
11 
16 
17 #include "CLHEP/Random/RandFlat.h"
18 
19 using namespace edm;
20 using namespace std;
21 
24 {
25 
26 
27  ParameterSet defpset ;
28  ParameterSet pgun_params =
29  pset.getParameter<ParameterSet>("PGunParameters") ;
30 
31  fMinPt = pgun_params.getParameter<double>("MinPt");
32  fMaxPt = pgun_params.getParameter<double>("MaxPt");
33 
34  fInConeIds = pgun_params.getParameter< vector<int> >("InConeID");
35  fMinDeltaR = pgun_params.getParameter<double>("MinDeltaR");
36  fMaxDeltaR = pgun_params.getParameter<double>("MaxDeltaR");
37  fMinMomRatio = pgun_params.getParameter<double>("MinMomRatio");
38  fMaxMomRatio = pgun_params.getParameter<double>("MaxMomRatio");
39 
40  fInConeMinEta = pgun_params.getParameter<double>("InConeMinEta");
41  fInConeMaxEta = pgun_params.getParameter<double>("InConeMaxEta");
42  fInConeMinPhi = pgun_params.getParameter<double>("InConeMinPhi");
43  fInConeMaxPhi = pgun_params.getParameter<double>("InConeMaxPhi");
44  fInConeMaxTry = pgun_params.getParameter<unsigned int>("InConeMaxTry");
45 
46  produces<HepMCProduct>("unsmeared");
47  produces<GenEventInfoProduct>();
48 }
49 
51 {
52  // no need to cleanup GenEvent memory - done in HepMCProduct
53 }
54 
56 {
58  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
59 
60  if ( fVerbosity > 0 )
61  {
62  cout << " MultiParticleInConeGunProducer : Begin New Event Generation" << endl ;
63  }
64  // event loop (well, another step in it...)
65 
66  // no need to clean up GenEvent memory - done in HepMCProduct
67  //
68 
69  // here re-create fEvt (memory)
70  //
71  fEvt = new HepMC::GenEvent() ;
72 
73  // now actualy, cook up the event from PDGTable and gun parameters
74  //
75  // 1st, primary vertex
76  //
77  //HepMC::GenVertex* Vtx = new HepMC::GenVertex(CLHEP::HepLorentzVector(0.,0.,0.));
78  HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
79 
80  // loop over particles
81  //
82  int barcode = 1 ;
83  for (unsigned int ip=0; ip<fPartIDs.size(); ++ip)
84  {
85 
86  double pt = CLHEP::RandFlat::shoot(engine, fMinPt, fMaxPt) ;
87  double eta = CLHEP::RandFlat::shoot(engine, fMinEta, fMaxEta) ;
88  double phi = CLHEP::RandFlat::shoot(engine, fMinPhi, fMaxPhi) ;
89  int PartID = fPartIDs[ip] ;
90  const HepPDT::ParticleData*
91  PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
92  double mass = PData->mass().value() ;
93  double theta = 2.*atan(exp(-eta)) ;
94  double mom = pt/sin(theta) ;
95  double px = pt*cos(phi) ;
96  double py = pt*sin(phi) ;
97  double pz = mom*cos(theta) ;
98  double energy2= mom*mom + mass*mass ;
99  double energy = sqrt(energy2) ;
100 
101  HepMC::FourVector p(px,py,pz,energy) ;
102  HepMC::GenParticle* Part =
103  new HepMC::GenParticle(p,PartID,1);
104  Part->suggest_barcode( barcode ) ;
105  barcode++ ;
106  Vtx->add_particle_out(Part);
107 
108  if ( fAddAntiParticle ){}
109 
110  // now add the particles in the cone
111  for (unsigned iPic=0; iPic!=fInConeIds.size();iPic++){
112  unsigned int nTry=0;
113  while(true){
114  //shoot flat Deltar
115  double dR = CLHEP::RandFlat::shoot(engine, fMinDeltaR, fMaxDeltaR);
116  //shoot flat eta/phi mixing
117  double alpha = CLHEP::RandFlat::shoot(engine, -3.14159265358979323846, 3.14159265358979323846);
118  double dEta = dR*cos(alpha);
119  double dPhi = dR*sin(alpha);
120 
121  /*
122  //shoot Energy of associated particle
123  double energyIc = CLHEP::RandFlat::shoot(engine, fMinEInCone, fMaxEInCone);
124  if (mom2Ic>0){ momIC = sqrt(mom2Ic);}
125  */
126  // get kinematics
127  double etaIc = eta+dEta;
128  double phiIc = phi+dPhi;
129  //put it back in -Pi:Pi if necessary. multiple time might be necessary if dR > 3
130  const unsigned int maxL=100;
131  unsigned int iL=0;
132  while(iL++<maxL){
133  if (phiIc > 3.14159265358979323846) phiIc-=2*3.14159265358979323846;
134  else if(phiIc <-3.14159265358979323846) phiIc+=2*3.14159265358979323846;
135 
136  if (abs(phiIc)<3.14159265358979323846) break;
137  }
138 
139 
140  //allow to skip it if you did not run out of possible drawing
141  if (nTry++<=fInConeMaxTry){
142  //draw another one if this one is not in acceptance
143  if (etaIc<fInConeMinEta || etaIc > fInConeMaxEta) continue;
144  if (phiIc<fInConeMinPhi || phiIc > fInConeMaxPhi) continue;
145  }
146  else{
147  if ( fVerbosity > 0 )
148  {
149  cout << " MultiParticleInConeGunProducer : could not produce a particle "
150  << fInConeIds[iPic]<<" in cone "
151  << fMinDeltaR<<" to "<<fMaxDeltaR<<" within the "<<fInConeMaxTry<<" allowed tries."<<endl;
152  }
153  /* edm::LogError("MultiParticleInConeGunProducer")<< " MultiParticleInConeGunProducer : could not produce a particle "<<
154  fInConeIds[iPic]<<" in cone "<<
155  fMinDeltaR<<" to "<<fMaxDeltaR<<" within the "<<fInConeMaxTry<<" allowed tries.";*/
156  }
157  int PartIDIc=fInConeIds[iPic];
158  const HepPDT::ParticleData*
159  PDataIc = fPDGTable->particle(HepPDT::ParticleID(abs(PartIDIc)));
160 
161  //shoot momentum ratio
162  double momR = CLHEP::RandFlat::shoot(engine, fMinMomRatio, fMaxMomRatio);
163  double massIc= PDataIc->mass().value() ;
164  double momIc = momR * mom;
165  double energyIc = sqrt(momIc*momIc + massIc*massIc);
166 
167  double thetaIc = 2.*atan(exp(-etaIc)) ;
168  double pxIc = momIc*sin(thetaIc)*cos(phiIc);
169  double pyIc = momIc*sin(thetaIc)*sin(phiIc);
170  double pzIc = momIc*cos(thetaIc);
171 
172  HepMC::FourVector pIc(pxIc,pyIc,pzIc,energyIc) ;
173  HepMC::GenParticle* PartIc = new HepMC::GenParticle(pIc, PartIDIc, 1);
174  PartIc->suggest_barcode( barcode ) ;
175  barcode++ ;
176  Vtx->add_particle_out(PartIc);
177  break;
178  }//try many times while not in acceptance
179  }//loop over the particle Ids in the cone
180  }
181 
182  fEvt->add_vertex(Vtx) ;
183  fEvt->set_event_number(e.id().event()) ;
184  fEvt->set_signal_process_id(20) ;
185 
186  if ( fVerbosity > 0 )
187  {
188  fEvt->print() ;
189  }
190 
191  auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
192  BProduct->addHepMCData( fEvt );
193  e.put(BProduct, "unsmeared");
194 
195  auto_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
196  e.put(genEventInfo);
197 
198  if ( fVerbosity > 0 )
199  {
200  cout << " MultiParticleInConeGunProducer : Event Generation Done " << endl;
201  }
202 }
203 
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
float alpha
Definition: AMPTWrapper.h:95
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
ESHandle< HepPDT::ParticleDataTable > fPDGTable
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
T sqrt(T t)
Definition: SSEVec.h:48
virtual void produce(Event &e, const EventSetup &es) override
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
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &) const =0
Use this engine in event methods.
edm::EventID id() const
Definition: EventBase.h:60
StreamID streamID() const
Definition: Event.h:79
tuple cout
Definition: gather_cfg.py:121