CMS 3D CMS Logo

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