CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ExpoRandomPGunProducer.cc
Go to the documentation of this file.
1 /*
2  * \author Jean-Roch Vlimant
3  * modified by S.Abdullin 04/02/2011
4  */
5 
6 #include <ostream>
7 
9 
12 
17 
18 #include "CLHEP/Random/RandFlat.h"
19 
20 using namespace edm;
21 using namespace std;
22 
24  ParameterSet defpset;
25  ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters");
26 
27  fMinP = pgun_params.getParameter<double>("MinP");
28  fMaxP = pgun_params.getParameter<double>("MaxP");
29 
30  produces<HepMCProduct>("unsmeared");
31  produces<GenEventInfoProduct>();
32 }
33 
35  // no need to cleanup GenEvent memory - done in HepMCProduct
36 }
37 
40  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
41 
42  if (fVerbosity > 0) {
43  std::cout << " ExpoRandomPGunProducer : Begin New Event Generation" << std::endl;
44  }
45  // event loop (well, another step in it...)
46 
47  // no need to clean up GenEvent memory - done in HepMCProduct
48  //
49 
50  // here re-create fEvt (memory)
51  //
52 
53  fEvt = new HepMC::GenEvent();
54 
55  // now actualy, cook up the event from PDGTable and gun parameters
56  //
57  // 1st, primary vertex
58  //
59  //HepMC::GenVertex* Vtx = new HepMC::GenVertex(CLHEP::HepLorentzVector(0.,0.,0.));
60  HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0., 0., 0.));
61 
62  // loop over particles
63  //
64  int barcode = 1;
65  for (unsigned int ip = 0; ip < fPartIDs.size(); ++ip) {
66  double pmom = CLHEP::RandFlat::shoot(engine, fMinP, fMaxP);
67  double y = (1. / fMinP) * CLHEP::RandFlat::shoot(engine, 0.0, 1.0);
68  double f = 1. / pmom;
69  bool accpt = (y < f);
70  //shoot until in the designated range
71  while ((pmom < fMinP || pmom > fMaxP) || !accpt) {
72  pmom = CLHEP::RandFlat::shoot(engine, fMinP, fMaxP);
73  y = (1. / fMinP) * CLHEP::RandFlat::shoot(engine, 0.0, 1.0);
74  f = 1. / pmom;
75  accpt = (y < f);
76  }
77 
78  double eta = CLHEP::RandFlat::shoot(engine, fMinEta, fMaxEta);
79  double phi = CLHEP::RandFlat::shoot(engine, fMinPhi, fMaxPhi);
80  int PartID = fPartIDs[ip];
81  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID)));
82  double mass = PData->mass().value();
83  double theta = 2. * atan(exp(-eta));
84  double mom = pmom;
85  double pt = mom * sin(theta);
86  double px = pt * cos(phi);
87  double py = pt * sin(phi);
88  double pz = mom * cos(theta);
89  double energy2 = mom * mom + mass * mass;
90  double energy = sqrt(energy2);
91  //CLHEP::Hep3Vector p(px,py,pz) ;
92  //HepMC::GenParticle* Part =
93  // new HepMC::GenParticle(CLHEP::HepLorentzVector(p,energy),PartID,1);
94  HepMC::FourVector p(px, py, pz, energy);
95  HepMC::GenParticle* Part = new HepMC::GenParticle(p, PartID, 1);
96  Part->suggest_barcode(barcode);
97  barcode++;
98  Vtx->add_particle_out(Part);
99 
100  if (fAddAntiParticle) {
101  //CLHEP::Hep3Vector ap(-px,-py,-pz) ;
102  HepMC::FourVector ap(-px, -py, -pz, energy);
103  int APartID = -PartID;
104  if (PartID == 22 || PartID == 23) {
105  APartID = PartID;
106  }
107  //HepMC::GenParticle* APart =
108  // new HepMC::GenParticle(CLHEP::HepLorentzVector(ap,energy),APartID,1);
109  HepMC::GenParticle* APart = new HepMC::GenParticle(ap, APartID, 1);
110  APart->suggest_barcode(barcode);
111  barcode++;
112  Vtx->add_particle_out(APart);
113  }
114  }
115 
116  fEvt->add_vertex(Vtx);
117  fEvt->set_event_number(e.id().event());
118  fEvt->set_signal_process_id(20);
119 
120  if (fVerbosity > 0) {
121  fEvt->print();
122  }
123 
124  unique_ptr<HepMCProduct> BProduct(new HepMCProduct());
125  BProduct->addHepMCData(fEvt);
126  e.put(std::move(BProduct), "unsmeared");
127 
128  unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
129  e.put(std::move(genEventInfo));
130 
131  if (fVerbosity > 0) {
132  // for testing purpose only
133  // fEvt->print() ; // prints empty info after it's made into edm::Event
134  std::cout << " FlatRandomPGunProducer : Event Generation Done " << std::endl;
135  }
136 }
EventNumber_t event() const
Definition: EventID.h:40
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
void produce(Event &e, const EventSetup &es) override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
ESHandle< HepPDT::ParticleDataTable > fPDGTable
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
def move
Definition: eostools.py:511
ExpoRandomPGunProducer(const ParameterSet &pset)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > fPartIDs
HepPDT::ParticleData ParticleData
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EventID id() const
Definition: EventBase.h:59
StreamID streamID() const
Definition: Event.h:98
tuple cout
Definition: gather_cfg.py:144