CMS 3D CMS Logo

CloseByParticleGunProducer.cc
Go to the documentation of this file.
1 #include <ostream>
2 #include <cmath>
3 
5 
8 
10 
16 
17 #include <CLHEP/Random/RandFlat.h>
18 #include <CLHEP/Units/SystemOfUnits.h>
19 #include <CLHEP/Units/GlobalPhysicalConstants.h>
20 #include <CLHEP/Random/RandFlat.h>
21 
22 using namespace edm;
23 using namespace std;
24 
26  : BaseFlatGunProducer(pset), m_fieldToken(esConsumes()) {
27  ParameterSet pgun_params = pset.getParameter<ParameterSet>("PGunParameters");
28  fControlledByEta = pgun_params.getParameter<bool>("ControlledByEta");
29  fControlledByREta = pgun_params.getParameter<bool>("ControlledByREta");
31  throw cms::Exception("CloseByParticleGunProducer")
32  << " Conflicting configuration, cannot have both ControlledByEta and ControlledByREta ";
33 
34  fVarMax = pgun_params.getParameter<double>("VarMax");
35  fVarMin = pgun_params.getParameter<double>("VarMin");
36  fMaxVarSpread = pgun_params.getParameter<bool>("MaxVarSpread");
37  fLogSpacedVar = pgun_params.getParameter<bool>("LogSpacedVar");
38  fFlatPtGeneration = pgun_params.getParameter<bool>("FlatPtGeneration");
39  if (fVarMin < 1 && !fFlatPtGeneration)
40  throw cms::Exception("CloseByParticleGunProducer")
41  << " Please choose a minimum energy greater than 1 GeV, otherwise time "
42  "information may be invalid or not reliable";
43  if (fVarMin < 0 && fLogSpacedVar)
44  throw cms::Exception("CloseByParticleGunProducer") << " Minimum energy must be greater than zero for log spacing";
45  else {
48  }
49 
51  fEtaMax = pgun_params.getParameter<double>("MaxEta");
52  fEtaMin = pgun_params.getParameter<double>("MinEta");
53  if (fEtaMax <= fEtaMin)
54  throw cms::Exception("CloseByParticleGunProducer") << " Please fix MinEta and MaxEta values in the configuration";
55  } else {
56  fRMax = pgun_params.getParameter<double>("RMax");
57  fRMin = pgun_params.getParameter<double>("RMin");
58  if (fRMax <= fRMin)
59  throw cms::Exception("CloseByParticleGunProducer") << " Please fix RMin and RMax values in the configuration";
60  }
61  if (!fControlledByREta) {
62  fZMax = pgun_params.getParameter<double>("ZMax");
63  fZMin = pgun_params.getParameter<double>("ZMin");
64 
65  if (fZMax <= fZMin)
66  throw cms::Exception("CloseByParticleGunProducer") << " Please fix ZMin and ZMax values in the configuration";
67  }
68  fDelta = pgun_params.getParameter<double>("Delta");
69  fPhiMin = pgun_params.getParameter<double>("MinPhi");
70  fPhiMax = pgun_params.getParameter<double>("MaxPhi");
71  fPointing = pgun_params.getParameter<bool>("Pointing");
72  fOverlapping = pgun_params.getParameter<bool>("Overlapping");
74  throw cms::Exception("CloseByParticleGunProducer")
75  << " Can't generate non pointing FlatPt samples; please disable FlatPt generation or generate pointing sample";
76  fRandomShoot = pgun_params.getParameter<bool>("RandomShoot");
77  fNParticles = pgun_params.getParameter<int>("NParticles");
78  fPartIDs = pgun_params.getParameter<vector<int>>("PartID");
79 
80  // set dt between particles
81  fUseDeltaT = pgun_params.getParameter<bool>("UseDeltaT");
82  fTMax = pgun_params.getParameter<double>("TMax");
83  fTMin = pgun_params.getParameter<double>("TMin");
84  if (fTMax <= fTMin)
85  throw cms::Exception("CloseByParticleGunProducer") << " Please fix TMin and TMax values in the configuration";
86  // set a fixed time offset for the particles
87  fOffsetFirst = pgun_params.getParameter<double>("OffsetFirst");
88 
89  produces<HepMCProduct>("unsmeared");
90  produces<GenEventInfoProduct>();
91 }
92 
94  // no need to cleanup GenEvent memory - done in HepMCProduct
95 }
96 
99  desc.add<bool>("AddAntiParticle", false);
100  {
102  psd0.add<bool>("ControlledByEta", false);
103  psd0.add<bool>("ControlledByREta", false);
104  psd0.add<double>("Delta", 10);
105  psd0.add<double>("VarMax", 200.0);
106  psd0.add<double>("VarMin", 25.0);
107  psd0.add<bool>("MaxVarSpread", false);
108  psd0.add<bool>("LogSpacedVar", false);
109  psd0.add<bool>("FlatPtGeneration", false);
110  psd0.add<double>("MaxEta", 2.7);
111  psd0.add<double>("MaxPhi", 3.14159265359);
112  psd0.add<double>("MinEta", 1.7);
113  psd0.add<double>("MinPhi", -3.14159265359);
114  psd0.add<int>("NParticles", 2);
115  psd0.add<bool>("Overlapping", false);
116  psd0.add<std::vector<int>>("PartID",
117  {
118  22,
119  });
120  psd0.add<bool>("Pointing", true);
121  psd0.add<double>("RMax", 120);
122  psd0.add<double>("RMin", 60);
123  psd0.add<bool>("RandomShoot", false);
124  psd0.add<double>("ZMax", 321);
125  psd0.add<double>("ZMin", 320);
126  psd0.add<bool>("UseDeltaT", false);
127  psd0.add<double>("TMin", 0.);
128  psd0.add<double>("TMax", 0.05);
129  psd0.add<double>("OffsetFirst", 0.);
130  desc.add<edm::ParameterSetDescription>("PGunParameters", psd0);
131  }
132  desc.addUntracked<int>("Verbosity", 0);
133  desc.addUntracked<unsigned int>("firstRun", 1);
134  desc.add<std::string>("psethack", "random particles in phi and r windows");
135  descriptions.add("CloseByParticleGunProducer", desc);
136 }
137 
140  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
141 
142  if (fVerbosity > 0) {
143  LogDebug("CloseByParticleGunProducer") << " CloseByParticleGunProducer : Begin New Event Generation" << endl;
144  }
145  fEvt = new HepMC::GenEvent();
146 
147  auto const& field = es.getData(m_fieldToken);
148 
149  int barcode = 1;
150  unsigned int numParticles = fRandomShoot ? CLHEP::RandFlat::shoot(engine, 1, fNParticles) : fNParticles;
151 
152  double phi = CLHEP::RandFlat::shoot(engine, fPhiMin, fPhiMax);
153  double fZ;
154  double fR, fEta;
155  double fT;
156 
157  if (!fControlledByREta) {
158  fZ = CLHEP::RandFlat::shoot(engine, fZMin, fZMax);
159 
160  if (!fControlledByEta) {
161  fR = CLHEP::RandFlat::shoot(engine, fRMin, fRMax);
162  fEta = asinh(fZ / fR);
163  } else {
164  fEta = CLHEP::RandFlat::shoot(engine, fEtaMin, fEtaMax);
165  fR = (fZ / sinh(fEta));
166  }
167  } else {
168  fR = CLHEP::RandFlat::shoot(engine, fRMin, fRMax);
169  fEta = CLHEP::RandFlat::shoot(engine, fEtaMin, fEtaMax);
170  fZ = sinh(fEta) / fR;
171  }
172 
173  if (fUseDeltaT) {
174  fT = CLHEP::RandFlat::shoot(engine, fTMin, fTMax);
175  } else {
176  fT = 0.;
177  }
178 
179  double tmpPhi = phi;
180  double tmpR = fR;
181 
182  // Loop over particles
183  for (unsigned int ip = 0; ip < numParticles; ++ip) {
184  if (fOverlapping) {
185  fR = CLHEP::RandFlat::shoot(engine, tmpR - fDelta, tmpR + fDelta);
186  phi = CLHEP::RandFlat::shoot(engine, tmpPhi - fDelta / fR, tmpPhi + fDelta / fR);
187  } else
188  phi += fDelta / fR;
189 
190  double fVar;
191  if (numParticles > 1 && fMaxVarSpread)
192  fVar = fVarMin + ip * (fVarMax - fVarMin) / (numParticles - 1);
193  else if (fLogSpacedVar) {
194  double fVar_log = CLHEP::RandFlat::shoot(engine, log_fVarMin, log_fVarMax);
195  fVar = std::exp(fVar_log);
196  }
197 
198  else
199  fVar = CLHEP::RandFlat::shoot(engine, fVarMin, fVarMax);
200 
201  int partIdx = CLHEP::RandFlat::shoot(engine, 0, fPartIDs.size());
202  int PartID = fPartIDs[partIdx];
203  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID)));
204  double mass = PData->mass().value();
205 
206  double mom, px, py, pz;
207  double energy;
208 
209  if (!fFlatPtGeneration) {
210  double mom2 = fVar * fVar - mass * mass;
211  mom = 0.;
212  if (mom2 > 0.) {
213  mom = sqrt(mom2);
214  }
215  px = 0.;
216  py = 0.;
217  pz = mom;
218  energy = fVar;
219  } else {
220  double theta = 2. * atan(exp(-fEta));
221  mom = fVar / sin(theta);
222  px = fVar * cos(phi);
223  py = fVar * sin(phi);
224  pz = mom * cos(theta);
225  double energy2 = mom * mom + mass * mass;
226  energy = sqrt(energy2);
227  }
228  // Compute Vertex Position
229  double x = fR * cos(phi);
230  double y = fR * sin(phi);
231 
232  HepMC::FourVector p(px, py, pz, energy);
233  // If we are requested to be pointing to (0,0,0), correct the momentum direction
234  if (fPointing) {
235  math::XYZVector direction(x, y, fZ);
236  math::XYZVector momentum = direction.unit() * mom;
237  p.setX(momentum.x());
238  p.setY(momentum.y());
239  p.setZ(momentum.z());
240  }
241 
242  // compute correct path assuming uniform magnetic field in CMS
243  double pathLength = 0.;
244  const double speed = p.pz() / p.e() * c_light / CLHEP::cm;
245  if (PData->charge()) {
246  // Radius [cm] = P[GeV/c] * 10^9 / (c[mm/ns] * 10^6 * q[C] * B[T]) * 100[cm/m]
247  const double radius = std::sqrt(p.px() * p.px() + p.py() * p.py()) * std::pow(10, 5) /
248  (c_light * field.inTesla({0.f, 0.f, 0.f}).z()); // cm
249  const double arc = 2 * asinf(std::sqrt(x * x + y * y) / (2 * radius)) * radius;
250  pathLength = std::sqrt(arc * arc + fZ * fZ);
251  } else {
252  pathLength = std::sqrt(x * x + y * y + fZ * fZ);
253  }
254 
255  // if not pointing time doesn't mean a lot, keep the old way
256  const double pathTime = fPointing ? (pathLength / speed) : (std::sqrt(x * x + y * y + fZ * fZ) / speed);
257  double timeOffset = fOffsetFirst + (pathTime + ip * fT) * CLHEP::ns * c_light;
258 
259  HepMC::GenVertex* Vtx =
260  new HepMC::GenVertex(HepMC::FourVector(x * CLHEP::cm, y * CLHEP::cm, fZ * CLHEP::cm, timeOffset));
261 
263  Part->suggest_barcode(barcode);
264  barcode++;
265 
266  Vtx->add_particle_out(Part);
267 
268  if (fVerbosity > 0) {
269  Vtx->print();
270  Part->print();
271  }
272  fEvt->add_vertex(Vtx);
273  }
274 
275  fEvt->set_event_number(e.id().event());
276  fEvt->set_signal_process_id(20);
277 
278  if (fVerbosity > 0) {
279  fEvt->print();
280  }
281 
282  unique_ptr<HepMCProduct> BProduct(new HepMCProduct());
283  BProduct->addHepMCData(fEvt);
284  e.put(std::move(BProduct), "unsmeared");
285 
286  unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
287  e.put(std::move(genEventInfo));
288 
289  if (fVerbosity > 0) {
290  LogDebug("CloseByParticleGunProducer") << " CloseByParticleGunProducer : Event Generation Done " << endl;
291  }
292 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
static void fillDescriptions(ConfigurationDescriptions &descriptions)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > m_fieldToken
ESHandle< HepPDT::ParticleDataTable > fPDGTable
T sqrt(T t)
Definition: SSEVec.h:23
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void produce(Event &e, const EventSetup &es) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
HepPDT::ParticleData ParticleData
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
float x
Geom::Theta< T > theta() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511
float fEta(int hwEta)
Definition: conversion.h:17
#define LogDebug(id)