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 
14 
15 // #include "FWCore/Utilities/interface/Exception.h"
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>();
47  produces<GenEventInfoProduct>();
48 }
49 
51 {
52  // no need to cleanup GenEvent memory - done in HepMCProduct
53 }
54 
56 {
57 
58  if ( fVerbosity > 0 )
59  {
60  cout << " MultiParticleInConeGunProducer : Begin New Event Generation" << endl ;
61  }
62  // event loop (well, another step in it...)
63 
64  // no need to clean up GenEvent memory - done in HepMCProduct
65  //
66 
67  // here re-create fEvt (memory)
68  //
69  fEvt = new HepMC::GenEvent() ;
70 
71  // now actualy, cook up the event from PDGTable and gun parameters
72  //
73  // 1st, primary vertex
74  //
75  //HepMC::GenVertex* Vtx = new HepMC::GenVertex(CLHEP::HepLorentzVector(0.,0.,0.));
76  HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
77 
78  // loop over particles
79  //
80  int barcode = 1 ;
81  for (unsigned int ip=0; ip<fPartIDs.size(); ++ip)
82  {
83 
84  double pt = fRandomGenerator->fire(fMinPt, fMaxPt) ;
85  double eta = fRandomGenerator->fire(fMinEta, fMaxEta) ;
86  double phi = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
87  int PartID = fPartIDs[ip] ;
88  const HepPDT::ParticleData*
89  PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
90  double mass = PData->mass().value() ;
91  double theta = 2.*atan(exp(-eta)) ;
92  double mom = pt/sin(theta) ;
93  double px = pt*cos(phi) ;
94  double py = pt*sin(phi) ;
95  double pz = mom*cos(theta) ;
96  double energy2= mom*mom + mass*mass ;
97  double energy = sqrt(energy2) ;
98 
99  HepMC::FourVector p(px,py,pz,energy) ;
100  HepMC::GenParticle* Part =
101  new HepMC::GenParticle(p,PartID,1);
102  Part->suggest_barcode( barcode ) ;
103  barcode++ ;
104  Vtx->add_particle_out(Part);
105 
106  if ( fAddAntiParticle ){}
107 
108  // now add the particles in the cone
109  for (unsigned iPic=0; iPic!=fInConeIds.size();iPic++){
110  unsigned int nTry=0;
111  while(true){
112  //shoot flat Deltar
113  double dR = fRandomGenerator->fire(fMinDeltaR, fMaxDeltaR);
114  //shoot flat eta/phi mixing
115  double alpha = fRandomGenerator->fire(-3.14159265358979323846, 3.14159265358979323846);
116  double dEta = dR*cos(alpha);
117  double dPhi = dR*sin(alpha);
118 
119  /*
120  //shoot Energy of associated particle
121  double energyIc = fRandomGenerator->fire(fMinEInCone, fMaxEInCone);
122  if (mom2Ic>0){ momIC = sqrt(mom2Ic);}
123  */
124  // get kinematics
125  double etaIc = eta+dEta;
126  double phiIc = phi+dPhi;
127  //put it back in -Pi:Pi if necessary. multiple time might be necessary if dR > 3
128  const unsigned int maxL=100;
129  unsigned int iL=0;
130  while(iL++<maxL){
131  if (phiIc > 3.14159265358979323846) phiIc-=2*3.14159265358979323846;
132  else if(phiIc <-3.14159265358979323846) phiIc+=2*3.14159265358979323846;
133 
134  if (abs(phiIc)<3.14159265358979323846) break;
135  }
136 
137 
138  //allow to skip it if you did not run out of possible drawing
139  if (nTry++<=fInConeMaxTry){
140  //draw another one if this one is not in acceptance
141  if (etaIc<fInConeMinEta || etaIc > fInConeMaxEta) continue;
142  if (phiIc<fInConeMinPhi || phiIc > fInConeMaxPhi) continue;
143  }
144  else{
145  if ( fVerbosity > 0 )
146  {
147  cout << " MultiParticleInConeGunProducer : could not produce a particle "
148  << fInConeIds[iPic]<<" in cone "
149  << fMinDeltaR<<" to "<<fMaxDeltaR<<" within the "<<fInConeMaxTry<<" allowed tries."<<endl;
150  }
151  /* edm::LogError("MultiParticleInConeGunProducer")<< " MultiParticleInConeGunProducer : could not produce a particle "<<
152  fInConeIds[iPic]<<" in cone "<<
153  fMinDeltaR<<" to "<<fMaxDeltaR<<" within the "<<fInConeMaxTry<<" allowed tries.";*/
154  }
155  int PartIDIc=fInConeIds[iPic];
156  const HepPDT::ParticleData*
157  PDataIc = fPDGTable->particle(HepPDT::ParticleID(abs(PartIDIc)));
158 
159  //shoot momentum ratio
160  double momR = fRandomGenerator->fire(fMinMomRatio, fMaxMomRatio);
161  double massIc= PDataIc->mass().value() ;
162  double momIc = momR * mom;
163  double energyIc = sqrt(momIc*momIc + massIc*massIc);
164 
165  double thetaIc = 2.*atan(exp(-etaIc)) ;
166  double pxIc = momIc*sin(thetaIc)*cos(phiIc);
167  double pyIc = momIc*sin(thetaIc)*sin(phiIc);
168  double pzIc = momIc*cos(thetaIc);
169 
170  HepMC::FourVector pIc(pxIc,pyIc,pzIc,energyIc) ;
171  HepMC::GenParticle* PartIc = new HepMC::GenParticle(pIc, PartIDIc, 1);
172  PartIc->suggest_barcode( barcode ) ;
173  barcode++ ;
174  Vtx->add_particle_out(PartIc);
175  break;
176  }//try many times while not in acceptance
177  }//loop over the particle Ids in the cone
178  }
179 
180  fEvt->add_vertex(Vtx) ;
181  fEvt->set_event_number(e.id().event()) ;
182  fEvt->set_signal_process_id(20) ;
183 
184  if ( fVerbosity > 0 )
185  {
186  fEvt->print() ;
187  }
188 
189  auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
190  BProduct->addHepMCData( fEvt );
191  e.put(BProduct);
192 
193  auto_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
194  e.put(genEventInfo);
195 
196  if ( fVerbosity > 0 )
197  {
198  cout << " MultiParticleInConeGunProducer : Event Generation Done " << endl;
199  }
200 }
201 
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
float alpha
Definition: AMPTWrapper.h:95
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
T eta() 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:116
T sqrt(T t)
Definition: SSEVec.h:48
virtual void produce(Event &e, const EventSetup &es) override
CLHEP::RandFlat * fRandomGenerator
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
edm::EventID id() const
Definition: EventBase.h:56
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10