CMS 3D CMS Logo

EmbeddingLHEProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: test/EmbeddingLHEProducer
4 // Class: EmbeddingLHEProducer
5 //
13 //
14 // Original Author: Stefan Wayand
15 // Created: Wed, 13 Jan 2016 08:15:01 GMT
16 //
17 //
18 
19 // system include files
20 #include <algorithm>
21 #include <iostream>
22 #include <iterator>
23 #include <fstream>
24 #include <string>
25 #include <memory>
26 #include "TLorentzVector.h"
27 
28 // user include files
31 
35 
37 
42 
48 
52 
57 #include "CLHEP/Random/RandExponential.h"
58 
59 //
60 // class declaration
61 //
62 
63 namespace CLHEP {
64  class HepRandomEngine;
65 }
66 
67 class EmbeddingLHEProducer : public edm::one::EDProducer<edm::BeginRunProducer, edm::EndRunProducer> {
68 public:
69  explicit EmbeddingLHEProducer(const edm::ParameterSet &);
70  ~EmbeddingLHEProducer() override;
71 
72  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
73 
74 private:
75  void beginJob() override;
76  void produce(edm::Event &, const edm::EventSetup &) override;
77  void endJob() override;
78 
79  void beginRunProduce(edm::Run &run, edm::EventSetup const &es) override;
80  void endRunProduce(edm::Run &, edm::EventSetup const &) override;
81 
82  void fill_lhe_from_mumu(TLorentzVector &positiveLepton,
83  TLorentzVector &negativeLepton,
84  lhef::HEPEUP &outlhe,
85  CLHEP::HepRandomEngine *engine);
86  void fill_lhe_with_particle(lhef::HEPEUP &outlhe, TLorentzVector &particle, int pdgid, double spin, double ctau);
87 
88  void transform_mumu_to_tautau(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton);
90  void assign_4vector(TLorentzVector &Lepton, const pat::Muon *muon, std::string FSRmode);
91  void mirror(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton);
92  void rotate180(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton);
93 
95 
100  const double tauMass_ = 1.77682;
101  const double elMass_ = 0.00051;
102  const int embeddingParticles[3]{11, 13, 15};
103 
104  std::ofstream file;
106 
107  // instead of reconstruted 4vectors of muons generated are taken to study FSR. Possible modes:
108  // afterFSR - muons without FSR photons taken into account
109  // beforeFSR - muons with FSR photons taken into account
111 };
112 
113 //
114 // constructors and destructor
115 //
117  //register your products
118  produces<LHEEventProduct>();
119  produces<LHERunInfoProduct, edm::Transition::BeginRun>();
120  produces<math::XYZTLorentzVectorD>("vertexPosition");
121 
122  muonsCollection_ = consumes<edm::View<pat::Muon>>(iConfig.getParameter<edm::InputTag>("src"));
123  vertexCollection_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
124  particleToEmbed_ = iConfig.getParameter<int>("particleToEmbed");
125  mirror_ = iConfig.getParameter<bool>("mirror");
126  rotate180_ = iConfig.getParameter<bool>("rotate180");
127  studyFSRmode_ = iConfig.getUntrackedParameter<std::string>("studyFSRmode", "");
128 
129  write_lheout = false;
130  std::string lhe_ouputfile = iConfig.getUntrackedParameter<std::string>("lhe_outputfilename", "");
131  if (!lhe_ouputfile.empty()) {
132  write_lheout = true;
133  file.open(lhe_ouputfile, std::fstream::out | std::fstream::trunc);
134  }
135 
136  //check if particle can be embedded
137  if (std::find(std::begin(embeddingParticles), std::end(embeddingParticles), particleToEmbed_) ==
138  std::end(embeddingParticles)) {
139  throw cms::Exception("Configuration") << "The given particle to embed is not in the list of allowed particles.";
140  }
141 
143  if (!rng.isAvailable()) {
144  throw cms::Exception("Configuration") << "The EmbeddingLHEProducer requires the RandomNumberGeneratorService\n"
145  "which is not present in the configuration file. \n"
146  "You must add the service\n"
147  "in the configuration file or remove the modules that require it.";
148  }
149 }
150 
152 
153 //
154 // member functions
155 //
156 
157 // ------------ method called to produce the data ------------
159  using namespace edm;
160 
162  CLHEP::HepRandomEngine *engine = &rng->getEngine(iEvent.streamID());
163 
165  iEvent.getByToken(muonsCollection_, muonHandle);
166  edm::View<pat::Muon> coll_muons = *muonHandle;
167 
168  Handle<std::vector<reco::Vertex>> coll_vertices;
169  iEvent.getByToken(vertexCollection_, coll_vertices);
170 
171  TLorentzVector positiveLepton, negativeLepton;
172  bool mu_plus_found = false;
173  bool mu_minus_found = false;
174  lhef::HEPEUP hepeup;
175  hepeup.IDPRUP = 0;
176  hepeup.XWGTUP = 1;
177  hepeup.SCALUP = -1;
178  hepeup.AQEDUP = -1;
179  hepeup.AQCDUP = -1;
180  // Assuming Pt-Order
181  for (edm::View<pat::Muon>::const_iterator muon = coll_muons.begin(); muon != coll_muons.end(); ++muon) {
182  if (muon->charge() == 1 && !mu_plus_found) {
183  assign_4vector(positiveLepton, &(*muon), studyFSRmode_);
184  mu_plus_found = true;
185  } else if (muon->charge() == -1 && !mu_minus_found) {
186  assign_4vector(negativeLepton, &(*muon), studyFSRmode_);
187  mu_minus_found = true;
188  } else if (mu_minus_found && mu_plus_found)
189  break;
190  }
191  mirror(positiveLepton, negativeLepton); // if no mirror, function does nothing.
192  rotate180(positiveLepton, negativeLepton); // if no rotate180, function does nothing
193  transform_mumu_to_tautau(positiveLepton, negativeLepton); // if MuonEmbedding, function does nothing.
194  fill_lhe_from_mumu(positiveLepton, negativeLepton, hepeup, engine);
195 
196  double originalXWGTUP_ = 1.;
197  std::unique_ptr<LHEEventProduct> product(new LHEEventProduct(hepeup, originalXWGTUP_));
198 
199  if (write_lheout)
200  std::copy(product->begin(), product->end(), std::ostream_iterator<std::string>(file));
201 
202  iEvent.put(std::move(product));
203  // Saving vertex position
204  std::unique_ptr<math::XYZTLorentzVectorD> vertex_position(
205  new math::XYZTLorentzVectorD(coll_vertices->at(0).x(), coll_vertices->at(0).y(), coll_vertices->at(0).z(), 0.0));
206  iEvent.put(std::move(vertex_position), "vertexPosition");
207 }
208 
209 // ------------ method called once each job just before starting event loop ------------
211 
212 // ------------ method called once each job just after ending the event loop ------------
214 
215 // ------------ method called when starting to processes a run ------------
216 
218  // fill HEPRUP common block and store in edm::Run
219  lhef::HEPRUP heprup;
220 
221  // set number of processes: 1 for Z to tau tau
222  heprup.resize(1);
223 
224  //Process independent information
225 
226  //beam particles ID (two protons)
227  //heprup.IDBMUP.first = 2212;
228  //heprup.IDBMUP.second = 2212;
229 
230  //beam particles energies (both 6.5 GeV)
231  //heprup.EBMUP.first = 6500.;
232  //heprup.EBMUP.second = 6500.;
233 
234  //take default pdf group for both beamparticles
235  //heprup.PDFGUP.first = -1;
236  //heprup.PDFGUP.second = -1;
237 
238  //take certan pdf set ID (same as in officially produced DYJets LHE files)
239  //heprup.PDFSUP.first = -1;
240  //heprup.PDFSUP.second = -1;
241 
242  //master switch for event weight iterpretation (same as in officially produced DYJets LHE files)
243  heprup.IDWTUP = 3;
244 
245  //Information for first process (Z to tau tau), for now only placeholder:
246  heprup.XSECUP[0] = 1.;
247  heprup.XERRUP[0] = 0;
248  heprup.XMAXUP[0] = 1;
249  heprup.LPRUP[0] = 1;
250 
251  std::unique_ptr<LHERunInfoProduct> runInfo(new LHERunInfoProduct(heprup));
252  runInfo->addHeader(give_slha());
253 
254  if (write_lheout)
255  std::copy(runInfo->begin(), runInfo->end(), std::ostream_iterator<std::string>(file));
256  run.put(std::move(runInfo));
257 }
258 
260  if (write_lheout) {
262  file.close();
263  }
264 }
265 
266 void EmbeddingLHEProducer::fill_lhe_from_mumu(TLorentzVector &positiveLepton,
267  TLorentzVector &negativeLepton,
268  lhef::HEPEUP &outlhe,
269  CLHEP::HepRandomEngine *engine) {
270  TLorentzVector Z = positiveLepton + negativeLepton;
271  int leptonPDGID = particleToEmbed_;
272 
273  // double tau_ctau = 0.00871100; //cm
274  double tau_ctau0 = 8.71100e-02; // mm (for Pythia)
275  double tau_ctau_p =
276  tau_ctau0 * CLHEP::RandExponential::shoot(engine); // return -std::log(HepRandom::getTheEngine()->flat());
277  // replaces tau = process[iNow].tau0() * rndmPtr->exp(); from pythia8212/src/ProcessContainer.cc which is not initialized for ProcessLevel:all = off mode (no beam particle mode)
278  double tau_ctau_n = tau_ctau0 * CLHEP::RandExponential::shoot(engine);
279  //std::cout<<"tau_ctau P: "<<tau_ctau_p<<" tau_ctau N: "<<tau_ctau_n<<std::endl;
280 
281  fill_lhe_with_particle(outlhe, Z, 23, 9.0, 0);
282  fill_lhe_with_particle(outlhe, positiveLepton, -leptonPDGID, 1.0, tau_ctau_p);
283  fill_lhe_with_particle(outlhe, negativeLepton, leptonPDGID, -1.0, tau_ctau_n);
284 
285  return;
286 }
287 
289  lhef::HEPEUP &outlhe, TLorentzVector &particle, int pdgid, double spin, double ctau) {
290  // Pay attention to different index conventions:
291  // 'particleindex' follows usual C++ index conventions starting at 0 for a list.
292  // 'motherindex' follows the LHE index conventions: 0 is for 'not defined', so the listing starts at 1.
293  // That means: LHE index 1 == C++ index 0.
294  int particleindex = outlhe.NUP;
295  outlhe.resize(outlhe.NUP + 1);
296 
297  outlhe.PUP[particleindex][0] = particle.Px();
298  outlhe.PUP[particleindex][1] = particle.Py();
299  outlhe.PUP[particleindex][2] = particle.Pz();
300  outlhe.PUP[particleindex][3] = particle.E();
301  outlhe.PUP[particleindex][4] = particle.M();
302  outlhe.IDUP[particleindex] = pdgid;
303  outlhe.SPINUP[particleindex] = spin;
304  outlhe.VTIMUP[particleindex] = ctau;
305 
306  outlhe.ICOLUP[particleindex].first = 0;
307  outlhe.ICOLUP[particleindex].second = 0;
308 
309  if (std::abs(pdgid) == 23) {
310  outlhe.MOTHUP[particleindex].first = 0; // No Mother
311  outlhe.MOTHUP[particleindex].second = 0;
312  outlhe.ISTUP[particleindex] = 2; // status
313  }
314 
315  if (std::find(std::begin(embeddingParticles), std::end(embeddingParticles), std::abs(pdgid)) !=
316  std::end(embeddingParticles)) {
317  outlhe.MOTHUP[particleindex].first = 1; // Mother is the Z (first partile)
318  outlhe.MOTHUP[particleindex].second = 1; // Mother is the Z (first partile)
319 
320  outlhe.ISTUP[particleindex] = 1; //status
321  }
322 
323  return;
324 }
325 
326 void EmbeddingLHEProducer::transform_mumu_to_tautau(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton) {
327  // No corrections applied for muon embedding
328  double lep_mass;
329  if (particleToEmbed_ == 11) {
330  lep_mass = elMass_;
331  } else if (particleToEmbed_ == 15) {
332  lep_mass = tauMass_;
333  } else {
334  return;
335  }
336 
337  TLorentzVector Z = positiveLepton + negativeLepton;
338 
339  TVector3 boost_from_Z_to_LAB = Z.BoostVector();
340  TVector3 boost_from_LAB_to_Z = -Z.BoostVector();
341 
342  // Boosting the two leptons to Z restframe, then both are back to back. This means, same 3-momentum squared
343  positiveLepton.Boost(boost_from_LAB_to_Z);
344  negativeLepton.Boost(boost_from_LAB_to_Z);
345 
346  // Energy of tau = 0.5*Z-mass
347  double lep_mass_squared = lep_mass * lep_mass;
348  double lep_energy_squared = 0.25 * Z.M2();
349  double lep_3momentum_squared = lep_energy_squared - lep_mass_squared;
350  if (lep_3momentum_squared < 0) {
351  edm::LogWarning("TauEmbedding") << "3-Momentum squared is negative";
352  return;
353  }
354 
355  //Computing scale, applying it on the 3-momenta and building new 4 momenta of the taus
356  double scale = std::sqrt(lep_3momentum_squared / positiveLepton.Vect().Mag2());
357  positiveLepton.SetPxPyPzE(scale * positiveLepton.Px(),
358  scale * positiveLepton.Py(),
359  scale * positiveLepton.Pz(),
360  std::sqrt(lep_energy_squared));
361  negativeLepton.SetPxPyPzE(scale * negativeLepton.Px(),
362  scale * negativeLepton.Py(),
363  scale * negativeLepton.Pz(),
364  std::sqrt(lep_energy_squared));
365 
366  //Boosting the new taus back to LAB frame
367  positiveLepton.Boost(boost_from_Z_to_LAB);
368  negativeLepton.Boost(boost_from_Z_to_LAB);
369 
370  return;
371 }
372 
373 void EmbeddingLHEProducer::assign_4vector(TLorentzVector &Lepton, const pat::Muon *muon, std::string FSRmode) {
374  if ("afterFSR" == FSRmode && muon->genParticle() != nullptr) {
375  const reco::GenParticle *afterFSRMuon = muon->genParticle();
376  Lepton.SetPxPyPzE(
377  afterFSRMuon->p4().px(), afterFSRMuon->p4().py(), afterFSRMuon->p4().pz(), afterFSRMuon->p4().e());
378  } else if ("beforeFSR" == FSRmode && muon->genParticle() != nullptr) {
379  const reco::Candidate *beforeFSRMuon = find_original_muon(muon->genParticle());
380  Lepton.SetPxPyPzE(
381  beforeFSRMuon->p4().px(), beforeFSRMuon->p4().py(), beforeFSRMuon->p4().pz(), beforeFSRMuon->p4().e());
382  } else {
383  Lepton.SetPxPyPzE(muon->p4().px(), muon->p4().py(), muon->p4().pz(), muon->p4().e());
384  }
385  return;
386 }
387 
389  if (muon->mother(0) == nullptr)
390  return muon;
391  if (muon->pdgId() == muon->mother(0)->pdgId())
392  return find_original_muon(muon->mother(0));
393  else
394  return muon;
395 }
396 
397 void EmbeddingLHEProducer::rotate180(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton) {
398  if (!rotate180_)
399  return;
400  edm::LogInfo("TauEmbedding") << "Applying 180<C2><B0> rotation";
401  // By construction, the 3-momenta of mu-, mu+ and Z are in one plane.
402  // That means, one vector for perpendicular projection can be used for both leptons.
403  TLorentzVector Z = positiveLepton + negativeLepton;
404 
405  edm::LogInfo("TauEmbedding") << "MuMinus before. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
406  << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
407 
408  TVector3 Z3 = Z.Vect();
409  TVector3 positiveLepton3 = positiveLepton.Vect();
410  TVector3 negativeLepton3 = negativeLepton.Vect();
411 
412  TVector3 p3_perp = positiveLepton3 - positiveLepton3.Dot(Z3) / Z3.Dot(Z3) * Z3;
413  p3_perp = p3_perp.Unit();
414 
415  positiveLepton3 -= 2 * positiveLepton3.Dot(p3_perp) * p3_perp;
416  negativeLepton3 -= 2 * negativeLepton3.Dot(p3_perp) * p3_perp;
417 
418  positiveLepton.SetVect(positiveLepton3);
419  negativeLepton.SetVect(negativeLepton3);
420 
421  edm::LogInfo("TauEmbedding") << "MuMinus after. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
422  << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
423 
424  return;
425 }
426 
427 void EmbeddingLHEProducer::mirror(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton) {
428  if (!mirror_)
429  return;
430  edm::LogInfo("TauEmbedding") << "Applying mirroring";
431  TLorentzVector Z = positiveLepton + negativeLepton;
432 
433  edm::LogInfo("TauEmbedding") << "MuMinus before. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
434  << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
435 
436  TVector3 Z3 = Z.Vect();
437  TVector3 positiveLepton3 = positiveLepton.Vect();
438  TVector3 negativeLepton3 = negativeLepton.Vect();
439 
440  TVector3 beam(0., 0., 1.);
441  TVector3 perpToZandBeam = Z3.Cross(beam).Unit();
442 
443  positiveLepton3 -= 2 * positiveLepton3.Dot(perpToZandBeam) * perpToZandBeam;
444  negativeLepton3 -= 2 * negativeLepton3.Dot(perpToZandBeam) * perpToZandBeam;
445 
446  positiveLepton.SetVect(positiveLepton3);
447  negativeLepton.SetVect(negativeLepton3);
448 
449  edm::LogInfo("TauEmbedding") << "MuMinus after. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
450  << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
451 
452  return;
453 }
454 
456  LHERunInfoProduct::Header slhah("slha");
457 
458  slhah.addLine("######################################################################\n");
459  slhah.addLine("## PARAM_CARD AUTOMATICALY GENERATED BY MG5 FOLLOWING UFO MODEL ####\n");
460  slhah.addLine("######################################################################\n");
461  slhah.addLine("## ##\n");
462  slhah.addLine("## Width set on Auto will be computed following the information ##\n");
463  slhah.addLine("## present in the decay.py files of the model. ##\n");
464  slhah.addLine("## See arXiv:1402.1178 for more details. ##\n");
465  slhah.addLine("## ##\n");
466  slhah.addLine("######################################################################\n");
467  slhah.addLine("\n");
468  slhah.addLine("###################################\n");
469  slhah.addLine("## INFORMATION FOR MASS\n");
470  slhah.addLine("###################################\n");
471  slhah.addLine("Block mass \n");
472  slhah.addLine(" 6 1.730000e+02 # MT \n");
473  slhah.addLine(" 15 1.777000e+00 # MTA \n");
474  slhah.addLine(" 23 9.118800e+01 # MZ \n");
475  slhah.addLine(" 25 1.250000e+02 # MH \n");
476  slhah.addLine("## Dependent parameters, given by model restrictions.\n");
477  slhah.addLine("## Those values should be edited following the \n");
478  slhah.addLine("## analytical expression. MG5 ignores those values \n");
479  slhah.addLine("## but they are important for interfacing the output of MG5\n");
480  slhah.addLine("## to external program such as Pythia.\n");
481  slhah.addLine(" 1 0.000000 # d : 0.0 \n");
482  slhah.addLine(" 2 0.000000 # u : 0.0 \n");
483  slhah.addLine(" 3 0.000000 # s : 0.0 \n");
484  slhah.addLine(" 4 0.000000 # c : 0.0 \n");
485  slhah.addLine(" 5 0.000000 # b : 0.0 \n");
486  slhah.addLine(" 11 0.000000 # e- : 0.0 \n");
487  slhah.addLine(" 12 0.000000 # ve : 0.0 \n");
488  slhah.addLine(" 13 0.000000 # mu- : 0.0 \n");
489  slhah.addLine(" 14 0.000000 # vm : 0.0 \n");
490  slhah.addLine(" 16 0.000000 # vt : 0.0 \n");
491  slhah.addLine(" 21 0.000000 # g : 0.0 \n");
492  slhah.addLine(" 22 0.000000 # a : 0.0 \n");
493  slhah.addLine(
494  " 24 80.419002 # w+ : cmath.sqrt(MZ__exp__2/2. + cmath.sqrt(MZ__exp__4/4. - "
495  "(aEW*cmath.pi*MZ__exp__2)/(Gf*sqrt__2))) \n");
496  slhah.addLine("\n");
497  slhah.addLine("###################################\n");
498  slhah.addLine("## INFORMATION FOR SMINPUTS\n");
499  slhah.addLine("###################################\n");
500  slhah.addLine("Block sminputs \n");
501  slhah.addLine(" 1 1.325070e+02 # aEWM1 \n");
502  slhah.addLine(" 2 1.166390e-05 # Gf \n");
503  slhah.addLine(" 3 1.180000e-01 # aS \n");
504  slhah.addLine("\n");
505  slhah.addLine("###################################\n");
506  slhah.addLine("## INFORMATION FOR WOLFENSTEIN\n");
507  slhah.addLine("###################################\n");
508  slhah.addLine("Block wolfenstein \n");
509  slhah.addLine(" 1 2.253000e-01 # lamWS \n");
510  slhah.addLine(" 2 8.080000e-01 # AWS \n");
511  slhah.addLine(" 3 1.320000e-01 # rhoWS \n");
512  slhah.addLine(" 4 3.410000e-01 # etaWS \n");
513  slhah.addLine("\n");
514  slhah.addLine("###################################\n");
515  slhah.addLine("## INFORMATION FOR YUKAWA\n");
516  slhah.addLine("###################################\n");
517  slhah.addLine("Block yukawa \n");
518  slhah.addLine(" 6 1.730000e+02 # ymt \n");
519  slhah.addLine(" 15 1.777000e+00 # ymtau \n");
520  slhah.addLine("\n");
521  slhah.addLine("###################################\n");
522  slhah.addLine("## INFORMATION FOR DECAY\n");
523  slhah.addLine("###################################\n");
524  slhah.addLine("DECAY 6 1.491500e+00 # WT \n");
525  slhah.addLine("DECAY 15 2.270000e-12 # WTau \n");
526  slhah.addLine("DECAY 23 2.441404e+00 # WZ \n");
527  slhah.addLine("DECAY 24 2.047600e+00 # WW \n");
528  slhah.addLine("DECAY 25 6.382339e-03 # WH \n");
529  slhah.addLine("## Dependent parameters, given by model restrictions.\n");
530  slhah.addLine("## Those values should be edited following the \n");
531  slhah.addLine("## analytical expression. MG5 ignores those values \n");
532  slhah.addLine("## but they are important for interfacing the output of MG5\n");
533  slhah.addLine("## to external program such as Pythia.\n");
534  slhah.addLine("DECAY 1 0.000000 # d : 0.0 \n");
535  slhah.addLine("DECAY 2 0.000000 # u : 0.0 \n");
536  slhah.addLine("DECAY 3 0.000000 # s : 0.0 \n");
537  slhah.addLine("DECAY 4 0.000000 # c : 0.0 \n");
538  slhah.addLine("DECAY 5 0.000000 # b : 0.0 \n");
539  slhah.addLine("DECAY 11 0.000000 # e- : 0.0 \n");
540  slhah.addLine("DECAY 12 0.000000 # ve : 0.0 \n");
541  slhah.addLine("DECAY 13 0.000000 # mu- : 0.0 \n");
542  slhah.addLine("DECAY 14 0.000000 # vm : 0.0 \n");
543  slhah.addLine("DECAY 16 0.000000 # vt : 0.0 \n");
544  slhah.addLine("DECAY 21 0.000000 # g : 0.0 \n");
545  slhah.addLine("DECAY 22 0.000000 # a : 0.0\n");
546 
547  return slhah;
548 }
549 
550 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
552  //The following says we do not know what parameters are allowed so do no validation
553  // Please change this to state exactly what you do use, even if it is no parameters
555  desc.setUnknown();
556  descriptions.addDefault(desc);
557 }
558 
559 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void resize(int nrup)
Definition: LesHouches.h:44
void mirror(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:240
std::vector< double > VTIMUP
Definition: LesHouches.h:252
void transform_mumu_to_tautau(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
const LorentzVector & p4() const final
four-momentum Lorentz vector
void fill_lhe_with_particle(lhef::HEPEUP &outlhe, TLorentzVector &particle, int pdgid, double spin, double ctau)
void beginRunProduce(edm::Run &run, edm::EventSetup const &es) override
T getUntrackedParameter(std::string const &, T const &) const
LHERunInfoProduct::Header give_slha()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
void resize(int nup)
Definition: LesHouches.h:161
void assign_4vector(TLorentzVector &Lepton, const pat::Muon *muon, std::string FSRmode)
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< FiveVector > PUP
Definition: LesHouches.h:246
std::vector< double > SPINUP
Definition: LesHouches.h:259
void addLine(const std::string &line)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< int > ISTUP
Definition: LesHouches.h:228
void endRunProduce(edm::Run &, edm::EventSetup const &) override
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:234
std::vector< int > IDUP
Definition: LesHouches.h:223
std::vector< double > XERRUP
Definition: LesHouches.h:118
Log< level::Info, false > LogInfo
std::vector< double > XMAXUP
Definition: LesHouches.h:123
double AQCDUP
Definition: LesHouches.h:218
void rotate180(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton)
edm::EDGetTokenT< edm::View< pat::Muon > > muonsCollection_
void produce(edm::Event &, const edm::EventSetup &) override
const reco::Candidate * find_original_muon(const reco::Candidate *muon)
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< reco::VertexCollection > vertexCollection_
static const std::string & endOfFile()
double AQEDUP
Definition: LesHouches.h:213
bool isAvailable() const
Definition: Service.h:40
void fill_lhe_from_mumu(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton, lhef::HEPEUP &outlhe, CLHEP::HepRandomEngine *engine)
Definition: Lepton.py:1
const_iterator begin() const
std::vector< double > XSECUP
Definition: LesHouches.h:112
Log< level::Warning, false > LogWarning
double XWGTUP
Definition: LesHouches.h:194
Analysis-level muon class.
Definition: Muon.h:51
const_iterator end() const
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
EmbeddingLHEProducer(const edm::ParameterSet &)
std::vector< int > LPRUP
Definition: LesHouches.h:128
double SCALUP
Definition: LesHouches.h:208
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector