CMS 3D CMS Logo

FlatRandomPtAndDxyGunProducer.cc
Go to the documentation of this file.
1 #include <ostream>
2 
4 
7 
12 
13 #include "CLHEP/Random/RandFlat.h"
14 
15 using namespace edm;
16 using namespace std;
17 
19  ParameterSet defpset;
20  ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters");
21 
22  fMinPt = pgun_params.getParameter<double>("MinPt");
23  fMaxPt = pgun_params.getParameter<double>("MaxPt");
24  dxyMin_ = pgun_params.getParameter<double>("dxyMin");
25  dxyMax_ = pgun_params.getParameter<double>("dxyMax");
26  lxyMax_ = pgun_params.getParameter<double>("LxyMax");
27  lzMax_ = pgun_params.getParameter<double>("LzMax");
28  ConeRadius_ = pgun_params.getParameter<double>("ConeRadius");
29  ConeH_ = pgun_params.getParameter<double>("ConeH");
30  DistanceToAPEX_ = pgun_params.getParameter<double>("DistanceToAPEX");
31 
32  produces<HepMCProduct>("unsmeared");
33  produces<GenEventInfoProduct>();
34 }
35 
37  // no need to cleanup GenEvent memory - done in HepMCProduct
38 }
39 
42  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
43 
44  if (fVerbosity > 0) {
45  cout << " FlatRandomPtAndDxyGunProducer : Begin New Event Generation" << endl;
46  }
47  // event loop (well, another step in it...)
48 
49  // no need to clean up GenEvent memory - done in HepMCProduct
50  //
51 
52  // here re-create fEvt (memory)
53  //
54  fEvt = new HepMC::GenEvent();
55 
56  // now actualy, cook up the event from PDGTable and gun parameters
57  int barcode = 1;
58  for (unsigned int ip = 0; ip < fPartIDs.size(); ++ip) {
59  double phi_vtx = 0;
60  double dxy = 0;
61  double pt = 0;
62  double eta = 0;
63  double px = 0;
64  double py = 0;
65  double pz = 0;
66  double vx = 0;
67  double vy = 0;
68  double vz = 0;
69  double lxy = 0;
70 
71  bool passLoop = false;
72  while (not passLoop) {
73  bool passLxy = false;
74  bool passLz = false;
75  phi_vtx = CLHEP::RandFlat::shoot(engine, fMinPhi, fMaxPhi);
76  dxy = CLHEP::RandFlat::shoot(engine, dxyMin_, dxyMax_);
77  float dxysign = CLHEP::RandFlat::shoot(engine, -1, 1);
78  if (dxysign < 0)
79  dxy = -dxy;
80  pt = CLHEP::RandFlat::shoot(engine, fMinPt, fMaxPt);
81  px = pt * cos(phi_vtx);
82  py = pt * sin(phi_vtx);
83  for (int i = 0; i < 10000; i++) {
84  vx = CLHEP::RandFlat::shoot(engine, -lxyMax_, lxyMax_);
85  vy = (pt * dxy + vx * py) / px;
86  lxy = sqrt(vx * vx + vy * vy);
87  if (lxy < abs(lxyMax_) and (vx * px + vy * py) > 0) {
88  passLxy = true;
89  break;
90  }
91  }
92  eta = CLHEP::RandFlat::shoot(engine, fMinEta, fMaxEta);
93  pz = pt * sinh(eta);
94  //vz = fabs(fRandomGaussGenerator->fire(0.0, LzWidth_/2.0));
95  float ConeTheta = ConeRadius_ / ConeH_;
96  for (int j = 0; j < 100; j++) {
97  vz = CLHEP::RandFlat::shoot(engine, 0.0, lzMax_); // this is abs(vz)
98  float v0 = vz - DistanceToAPEX_;
99  if (v0 <= 0 or lxy * lxy / (ConeTheta * ConeTheta) > v0 * v0) {
100  passLz = true;
101  break;
102  }
103  }
104  if (pz < 0)
105  vz = -vz;
106  passLoop = (passLxy and passLz);
107 
108  if (passLoop)
109  break;
110  }
111 
112  HepMC::GenVertex* Vtx1 = new HepMC::GenVertex(HepMC::FourVector(vx, vy, vz));
113 
114  int PartID = fPartIDs[ip];
115  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID)));
116  double mass = PData->mass().value();
117  double energy2 = px * px + py * py + pz * pz + mass * mass;
118  double energy = sqrt(energy2);
119  HepMC::FourVector p(px, py, pz, energy);
120  HepMC::GenParticle* Part = new HepMC::GenParticle(p, PartID, 1);
121  Part->suggest_barcode(barcode);
122  barcode++;
123  Vtx1->add_particle_out(Part);
124  fEvt->add_vertex(Vtx1);
125 
126  if (fAddAntiParticle) {
127  HepMC::GenVertex* Vtx2 = new HepMC::GenVertex(HepMC::FourVector(-vx, -vy, -vz));
128  HepMC::FourVector ap(-px, -py, -pz, energy);
129  int APartID = -PartID;
130  if (PartID == 22 || PartID == 23) {
131  APartID = PartID;
132  }
133  HepMC::GenParticle* APart = new HepMC::GenParticle(ap, APartID, 1);
134  APart->suggest_barcode(barcode);
135  barcode++;
136  Vtx2->add_particle_out(APart);
137  fEvt->add_vertex(Vtx2);
138  }
139  }
140  fEvt->set_event_number(e.id().event());
141  fEvt->set_signal_process_id(20);
142 
143  if (fVerbosity > 0) {
144  fEvt->print();
145  }
146 
147  unique_ptr<HepMCProduct> BProduct(new HepMCProduct());
148  BProduct->addHepMCData(fEvt);
149  e.put(std::move(BProduct), "unsmeared");
150 
151  unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
152  e.put(std::move(genEventInfo));
153 
154  if (fVerbosity > 0) {
155  cout << " FlatRandomPtAndDxyGunProducer : End New Event Generation" << endl;
156  fEvt->print();
157  }
158 }
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:40
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
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:19
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
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
FlatRandomPtAndDxyGunProducer(const ParameterSet &pset)
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
StreamID streamID() const
Definition: Event.h:96
def move(src, dest)
Definition: eostools.py:511