Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:25:59

0001 // -*- C++ -*-
0002 //
0003 // Package:    test/EmbeddingLHEProducer
0004 // Class:      EmbeddingLHEProducer
0005 //
0006 /**\class EmbeddingLHEProducer EmbeddingLHEProducer.cc test/EmbeddingLHEProducer/plugins/EmbeddingLHEProducer.cc
0007 
0008  Description: [one line class summary]
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  Stefan Wayand
0015 //         Created:  Wed, 13 Jan 2016 08:15:01 GMT
0016 //
0017 //
0018 
0019 // system include files
0020 #include <algorithm>
0021 #include <iostream>
0022 #include <iterator>
0023 #include <fstream>
0024 #include <string>
0025 #include <memory>
0026 #include "TLorentzVector.h"
0027 
0028 // user include files
0029 #include "FWCore/Framework/interface/Frameworkfwd.h"
0030 #include "FWCore/Framework/interface/one/EDProducer.h"
0031 
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/Run.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035 
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 
0038 #include "DataFormats/VertexReco/interface/Vertex.h"
0039 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0040 #include "DataFormats/PatCandidates/interface/Muon.h"
0041 #include "DataFormats/Math/interface/LorentzVector.h"
0042 
0043 #include "SimDataFormats/GeneratorProducts/interface/LesHouches.h"
0044 #include "SimDataFormats/GeneratorProducts/interface/LHECommonBlocks.h"
0045 #include "SimDataFormats/GeneratorProducts/interface/LHERunInfoProduct.h"
0046 #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
0047 #include "SimDataFormats/GeneratorProducts/interface/LHEXMLStringProduct.h"
0048 
0049 #include "GeneratorInterface/LHEInterface/interface/LHERunInfo.h"
0050 #include "GeneratorInterface/LHEInterface/interface/LHEEvent.h"
0051 #include "GeneratorInterface/LHEInterface/interface/LHEReader.h"
0052 
0053 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0054 #include "FWCore/ServiceRegistry/interface/Service.h"
0055 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
0056 #include "FWCore/Utilities/interface/StreamID.h"
0057 #include "CLHEP/Random/RandExponential.h"
0058 
0059 //
0060 // class declaration
0061 //
0062 
0063 namespace CLHEP {
0064   class HepRandomEngine;
0065 }
0066 
0067 class EmbeddingLHEProducer : public edm::one::EDProducer<edm::BeginRunProducer, edm::EndRunProducer> {
0068 public:
0069   explicit EmbeddingLHEProducer(const edm::ParameterSet &);
0070   ~EmbeddingLHEProducer() override;
0071 
0072   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0073 
0074 private:
0075   void beginJob() override;
0076   void produce(edm::Event &, const edm::EventSetup &) override;
0077   void endJob() override;
0078 
0079   void beginRunProduce(edm::Run &run, edm::EventSetup const &es) override;
0080   void endRunProduce(edm::Run &, edm::EventSetup const &) override;
0081 
0082   void fill_lhe_from_mumu(TLorentzVector &positiveLepton,
0083                           TLorentzVector &negativeLepton,
0084                           lhef::HEPEUP &outlhe,
0085                           CLHEP::HepRandomEngine *engine);
0086   void fill_lhe_with_particle(lhef::HEPEUP &outlhe, TLorentzVector &particle, int pdgid, double spin, double ctau);
0087 
0088   void transform_mumu_to_tautau(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton);
0089   const reco::Candidate *find_original_muon(const reco::Candidate *muon);
0090   void assign_4vector(TLorentzVector &Lepton, const pat::Muon *muon, std::string FSRmode);
0091   void mirror(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton);
0092   void rotate180(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton);
0093 
0094   LHERunInfoProduct::Header give_slha();
0095 
0096   edm::EDGetTokenT<edm::View<pat::Muon>> muonsCollection_;
0097   edm::EDGetTokenT<reco::VertexCollection> vertexCollection_;
0098   int particleToEmbed_;
0099   bool mirror_, rotate180_;
0100   const double tauMass_ = 1.77682;
0101   const double elMass_ = 0.00051;
0102   const int embeddingParticles[3]{11, 13, 15};
0103 
0104   std::ofstream file;
0105   bool write_lheout;
0106 
0107   // instead of reconstruted 4vectors of muons generated are taken to study FSR. Possible modes:
0108   // afterFSR - muons without FSR photons taken into account
0109   // beforeFSR - muons with FSR photons taken into account
0110   std::string studyFSRmode_;
0111 };
0112 
0113 //
0114 // constructors and destructor
0115 //
0116 EmbeddingLHEProducer::EmbeddingLHEProducer(const edm::ParameterSet &iConfig) {
0117   //register your products
0118   produces<LHEEventProduct>();
0119   produces<LHERunInfoProduct, edm::Transition::BeginRun>();
0120   produces<math::XYZTLorentzVectorD>("vertexPosition");
0121 
0122   muonsCollection_ = consumes<edm::View<pat::Muon>>(iConfig.getParameter<edm::InputTag>("src"));
0123   vertexCollection_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
0124   particleToEmbed_ = iConfig.getParameter<int>("particleToEmbed");
0125   mirror_ = iConfig.getParameter<bool>("mirror");
0126   rotate180_ = iConfig.getParameter<bool>("rotate180");
0127   studyFSRmode_ = iConfig.getUntrackedParameter<std::string>("studyFSRmode", "");
0128 
0129   write_lheout = false;
0130   std::string lhe_ouputfile = iConfig.getUntrackedParameter<std::string>("lhe_outputfilename", "");
0131   if (!lhe_ouputfile.empty()) {
0132     write_lheout = true;
0133     file.open(lhe_ouputfile, std::fstream::out | std::fstream::trunc);
0134   }
0135 
0136   //check if particle can be embedded
0137   if (std::find(std::begin(embeddingParticles), std::end(embeddingParticles), particleToEmbed_) ==
0138       std::end(embeddingParticles)) {
0139     throw cms::Exception("Configuration") << "The given particle to embed is not in the list of allowed particles.";
0140   }
0141 
0142   edm::Service<edm::RandomNumberGenerator> rng;
0143   if (!rng.isAvailable()) {
0144     throw cms::Exception("Configuration") << "The EmbeddingLHEProducer requires the RandomNumberGeneratorService\n"
0145                                              "which is not present in the configuration file. \n"
0146                                              "You must add the service\n"
0147                                              "in the configuration file or remove the modules that require it.";
0148   }
0149 }
0150 
0151 EmbeddingLHEProducer::~EmbeddingLHEProducer() {}
0152 
0153 //
0154 // member functions
0155 //
0156 
0157 // ------------ method called to produce the data  ------------
0158 void EmbeddingLHEProducer::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0159   using namespace edm;
0160 
0161   edm::Service<edm::RandomNumberGenerator> rng;
0162   CLHEP::HepRandomEngine *engine = &rng->getEngine(iEvent.streamID());
0163 
0164   edm::Handle<edm::View<pat::Muon>> muonHandle;
0165   iEvent.getByToken(muonsCollection_, muonHandle);
0166   edm::View<pat::Muon> coll_muons = *muonHandle;
0167 
0168   Handle<std::vector<reco::Vertex>> coll_vertices;
0169   iEvent.getByToken(vertexCollection_, coll_vertices);
0170 
0171   TLorentzVector positiveLepton, negativeLepton;
0172   bool mu_plus_found = false;
0173   bool mu_minus_found = false;
0174   lhef::HEPEUP hepeup;
0175   hepeup.IDPRUP = 0;
0176   hepeup.XWGTUP = 1;
0177   hepeup.SCALUP = -1;
0178   hepeup.AQEDUP = -1;
0179   hepeup.AQCDUP = -1;
0180   // Assuming Pt-Order
0181   for (edm::View<pat::Muon>::const_iterator muon = coll_muons.begin(); muon != coll_muons.end(); ++muon) {
0182     if (muon->charge() == 1 && !mu_plus_found) {
0183       assign_4vector(positiveLepton, &(*muon), studyFSRmode_);
0184       mu_plus_found = true;
0185     } else if (muon->charge() == -1 && !mu_minus_found) {
0186       assign_4vector(negativeLepton, &(*muon), studyFSRmode_);
0187       mu_minus_found = true;
0188     } else if (mu_minus_found && mu_plus_found)
0189       break;
0190   }
0191   mirror(positiveLepton, negativeLepton);                    // if no mirror, function does nothing.
0192   rotate180(positiveLepton, negativeLepton);                 // if no rotate180, function does nothing
0193   transform_mumu_to_tautau(positiveLepton, negativeLepton);  // if MuonEmbedding, function does nothing.
0194   fill_lhe_from_mumu(positiveLepton, negativeLepton, hepeup, engine);
0195 
0196   double originalXWGTUP_ = 1.;
0197   std::unique_ptr<LHEEventProduct> product(new LHEEventProduct(hepeup, originalXWGTUP_));
0198 
0199   if (write_lheout)
0200     std::copy(product->begin(), product->end(), std::ostream_iterator<std::string>(file));
0201 
0202   iEvent.put(std::move(product));
0203   // Saving vertex position
0204   std::unique_ptr<math::XYZTLorentzVectorD> vertex_position(
0205       new math::XYZTLorentzVectorD(coll_vertices->at(0).x(), coll_vertices->at(0).y(), coll_vertices->at(0).z(), 0.0));
0206   iEvent.put(std::move(vertex_position), "vertexPosition");
0207 }
0208 
0209 // ------------ method called once each job just before starting event loop  ------------
0210 void EmbeddingLHEProducer::beginJob() {}
0211 
0212 // ------------ method called once each job just after ending the event loop  ------------
0213 void EmbeddingLHEProducer::endJob() {}
0214 
0215 // ------------ method called when starting to processes a run  ------------
0216 
0217 void EmbeddingLHEProducer::beginRunProduce(edm::Run &run, edm::EventSetup const &) {
0218   // fill HEPRUP common block and store in edm::Run
0219   lhef::HEPRUP heprup;
0220 
0221   // set number of processes: 1 for Z to tau tau
0222   heprup.resize(1);
0223 
0224   //Process independent information
0225 
0226   //beam particles ID (two protons)
0227   //heprup.IDBMUP.first = 2212;
0228   //heprup.IDBMUP.second = 2212;
0229 
0230   //beam particles energies (both 6.5 GeV)
0231   //heprup.EBMUP.first = 6500.;
0232   //heprup.EBMUP.second = 6500.;
0233 
0234   //take default pdf group for both beamparticles
0235   //heprup.PDFGUP.first = -1;
0236   //heprup.PDFGUP.second = -1;
0237 
0238   //take certan pdf set ID (same as in officially produced DYJets LHE files)
0239   //heprup.PDFSUP.first = -1;
0240   //heprup.PDFSUP.second = -1;
0241 
0242   //master switch for event weight iterpretation (same as in officially produced DYJets LHE files)
0243   heprup.IDWTUP = 3;
0244 
0245   //Information for first process (Z to tau tau), for now only placeholder:
0246   heprup.XSECUP[0] = 1.;
0247   heprup.XERRUP[0] = 0;
0248   heprup.XMAXUP[0] = 1;
0249   heprup.LPRUP[0] = 1;
0250 
0251   std::unique_ptr<LHERunInfoProduct> runInfo(new LHERunInfoProduct(heprup));
0252   runInfo->addHeader(give_slha());
0253 
0254   if (write_lheout)
0255     std::copy(runInfo->begin(), runInfo->end(), std::ostream_iterator<std::string>(file));
0256   run.put(std::move(runInfo));
0257 }
0258 
0259 void EmbeddingLHEProducer::endRunProduce(edm::Run &run, edm::EventSetup const &es) {
0260   if (write_lheout) {
0261     file << LHERunInfoProduct::endOfFile();
0262     file.close();
0263   }
0264 }
0265 
0266 void EmbeddingLHEProducer::fill_lhe_from_mumu(TLorentzVector &positiveLepton,
0267                                               TLorentzVector &negativeLepton,
0268                                               lhef::HEPEUP &outlhe,
0269                                               CLHEP::HepRandomEngine *engine) {
0270   TLorentzVector Z = positiveLepton + negativeLepton;
0271   int leptonPDGID = particleToEmbed_;
0272 
0273   // double tau_ctau = 0.00871100; //cm
0274   double tau_ctau0 = 8.71100e-02;  // mm (for Pythia)
0275   double tau_ctau_p =
0276       tau_ctau0 * CLHEP::RandExponential::shoot(engine);  // return -std::log(HepRandom::getTheEngine()->flat());
0277   // 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)
0278   double tau_ctau_n = tau_ctau0 * CLHEP::RandExponential::shoot(engine);
0279   //std::cout<<"tau_ctau P: "<<tau_ctau_p<<" tau_ctau N:  "<<tau_ctau_n<<std::endl;
0280 
0281   fill_lhe_with_particle(outlhe, Z, 23, 9.0, 0);
0282   fill_lhe_with_particle(outlhe, positiveLepton, -leptonPDGID, 1.0, tau_ctau_p);
0283   fill_lhe_with_particle(outlhe, negativeLepton, leptonPDGID, -1.0, tau_ctau_n);
0284 
0285   return;
0286 }
0287 
0288 void EmbeddingLHEProducer::fill_lhe_with_particle(
0289     lhef::HEPEUP &outlhe, TLorentzVector &particle, int pdgid, double spin, double ctau) {
0290   // Pay attention to different index conventions:
0291   // 'particleindex' follows usual C++ index conventions starting at 0 for a list.
0292   // 'motherindex' follows the LHE index conventions: 0 is for 'not defined', so the listing starts at 1.
0293   // That means: LHE index 1 == C++ index 0.
0294   int particleindex = outlhe.NUP;
0295   outlhe.resize(outlhe.NUP + 1);
0296 
0297   outlhe.PUP[particleindex][0] = particle.Px();
0298   outlhe.PUP[particleindex][1] = particle.Py();
0299   outlhe.PUP[particleindex][2] = particle.Pz();
0300   outlhe.PUP[particleindex][3] = particle.E();
0301   outlhe.PUP[particleindex][4] = particle.M();
0302   outlhe.IDUP[particleindex] = pdgid;
0303   outlhe.SPINUP[particleindex] = spin;
0304   outlhe.VTIMUP[particleindex] = ctau;
0305 
0306   outlhe.ICOLUP[particleindex].first = 0;
0307   outlhe.ICOLUP[particleindex].second = 0;
0308 
0309   if (std::abs(pdgid) == 23) {
0310     outlhe.MOTHUP[particleindex].first = 0;  // No Mother
0311     outlhe.MOTHUP[particleindex].second = 0;
0312     outlhe.ISTUP[particleindex] = 2;  // status
0313   }
0314 
0315   if (std::find(std::begin(embeddingParticles), std::end(embeddingParticles), std::abs(pdgid)) !=
0316       std::end(embeddingParticles)) {
0317     outlhe.MOTHUP[particleindex].first = 1;   // Mother is the Z (first partile)
0318     outlhe.MOTHUP[particleindex].second = 1;  // Mother is the Z (first partile)
0319 
0320     outlhe.ISTUP[particleindex] = 1;  //status
0321   }
0322 
0323   return;
0324 }
0325 
0326 void EmbeddingLHEProducer::transform_mumu_to_tautau(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton) {
0327   // No corrections applied for muon embedding
0328   double lep_mass;
0329   if (particleToEmbed_ == 11) {
0330     lep_mass = elMass_;
0331   } else if (particleToEmbed_ == 15) {
0332     lep_mass = tauMass_;
0333   } else {
0334     return;
0335   }
0336 
0337   TLorentzVector Z = positiveLepton + negativeLepton;
0338 
0339   TVector3 boost_from_Z_to_LAB = Z.BoostVector();
0340   TVector3 boost_from_LAB_to_Z = -Z.BoostVector();
0341 
0342   // Boosting the two leptons to Z restframe, then both are back to back. This means, same 3-momentum squared
0343   positiveLepton.Boost(boost_from_LAB_to_Z);
0344   negativeLepton.Boost(boost_from_LAB_to_Z);
0345 
0346   // Energy of tau = 0.5*Z-mass
0347   double lep_mass_squared = lep_mass * lep_mass;
0348   double lep_energy_squared = 0.25 * Z.M2();
0349   double lep_3momentum_squared = lep_energy_squared - lep_mass_squared;
0350   if (lep_3momentum_squared < 0) {
0351     edm::LogWarning("TauEmbedding") << "3-Momentum squared is negative";
0352     return;
0353   }
0354 
0355   //Computing scale, applying it on the 3-momenta and building new 4 momenta of the taus
0356   double scale = std::sqrt(lep_3momentum_squared / positiveLepton.Vect().Mag2());
0357   positiveLepton.SetPxPyPzE(scale * positiveLepton.Px(),
0358                             scale * positiveLepton.Py(),
0359                             scale * positiveLepton.Pz(),
0360                             std::sqrt(lep_energy_squared));
0361   negativeLepton.SetPxPyPzE(scale * negativeLepton.Px(),
0362                             scale * negativeLepton.Py(),
0363                             scale * negativeLepton.Pz(),
0364                             std::sqrt(lep_energy_squared));
0365 
0366   //Boosting the new taus back to LAB frame
0367   positiveLepton.Boost(boost_from_Z_to_LAB);
0368   negativeLepton.Boost(boost_from_Z_to_LAB);
0369 
0370   return;
0371 }
0372 
0373 void EmbeddingLHEProducer::assign_4vector(TLorentzVector &Lepton, const pat::Muon *muon, std::string FSRmode) {
0374   if ("afterFSR" == FSRmode && muon->genParticle() != nullptr) {
0375     const reco::GenParticle *afterFSRMuon = muon->genParticle();
0376     Lepton.SetPxPyPzE(
0377         afterFSRMuon->p4().px(), afterFSRMuon->p4().py(), afterFSRMuon->p4().pz(), afterFSRMuon->p4().e());
0378   } else if ("beforeFSR" == FSRmode && muon->genParticle() != nullptr) {
0379     const reco::Candidate *beforeFSRMuon = find_original_muon(muon->genParticle());
0380     Lepton.SetPxPyPzE(
0381         beforeFSRMuon->p4().px(), beforeFSRMuon->p4().py(), beforeFSRMuon->p4().pz(), beforeFSRMuon->p4().e());
0382   } else {
0383     Lepton.SetPxPyPzE(muon->p4().px(), muon->p4().py(), muon->p4().pz(), muon->p4().e());
0384   }
0385   return;
0386 }
0387 
0388 const reco::Candidate *EmbeddingLHEProducer::find_original_muon(const reco::Candidate *muon) {
0389   if (muon->mother(0) == nullptr)
0390     return muon;
0391   if (muon->pdgId() == muon->mother(0)->pdgId())
0392     return find_original_muon(muon->mother(0));
0393   else
0394     return muon;
0395 }
0396 
0397 void EmbeddingLHEProducer::rotate180(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton) {
0398   if (!rotate180_)
0399     return;
0400   edm::LogInfo("TauEmbedding") << "Applying 180<C2><B0> rotation";
0401   // By construction, the 3-momenta of mu-, mu+ and Z are in one plane.
0402   // That means, one vector for perpendicular projection can be used for both leptons.
0403   TLorentzVector Z = positiveLepton + negativeLepton;
0404 
0405   edm::LogInfo("TauEmbedding") << "MuMinus before. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
0406                                << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
0407 
0408   TVector3 Z3 = Z.Vect();
0409   TVector3 positiveLepton3 = positiveLepton.Vect();
0410   TVector3 negativeLepton3 = negativeLepton.Vect();
0411 
0412   TVector3 p3_perp = positiveLepton3 - positiveLepton3.Dot(Z3) / Z3.Dot(Z3) * Z3;
0413   p3_perp = p3_perp.Unit();
0414 
0415   positiveLepton3 -= 2 * positiveLepton3.Dot(p3_perp) * p3_perp;
0416   negativeLepton3 -= 2 * negativeLepton3.Dot(p3_perp) * p3_perp;
0417 
0418   positiveLepton.SetVect(positiveLepton3);
0419   negativeLepton.SetVect(negativeLepton3);
0420 
0421   edm::LogInfo("TauEmbedding") << "MuMinus after. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
0422                                << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
0423 
0424   return;
0425 }
0426 
0427 void EmbeddingLHEProducer::mirror(TLorentzVector &positiveLepton, TLorentzVector &negativeLepton) {
0428   if (!mirror_)
0429     return;
0430   edm::LogInfo("TauEmbedding") << "Applying mirroring";
0431   TLorentzVector Z = positiveLepton + negativeLepton;
0432 
0433   edm::LogInfo("TauEmbedding") << "MuMinus before. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
0434                                << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
0435 
0436   TVector3 Z3 = Z.Vect();
0437   TVector3 positiveLepton3 = positiveLepton.Vect();
0438   TVector3 negativeLepton3 = negativeLepton.Vect();
0439 
0440   TVector3 beam(0., 0., 1.);
0441   TVector3 perpToZandBeam = Z3.Cross(beam).Unit();
0442 
0443   positiveLepton3 -= 2 * positiveLepton3.Dot(perpToZandBeam) * perpToZandBeam;
0444   negativeLepton3 -= 2 * negativeLepton3.Dot(perpToZandBeam) * perpToZandBeam;
0445 
0446   positiveLepton.SetVect(positiveLepton3);
0447   negativeLepton.SetVect(negativeLepton3);
0448 
0449   edm::LogInfo("TauEmbedding") << "MuMinus after. Pt: " << negativeLepton.Pt() << " Eta: " << negativeLepton.Eta()
0450                                << " Phi: " << negativeLepton.Phi() << " Mass: " << negativeLepton.M();
0451 
0452   return;
0453 }
0454 
0455 LHERunInfoProduct::Header EmbeddingLHEProducer::give_slha() {
0456   LHERunInfoProduct::Header slhah("slha");
0457 
0458   slhah.addLine("######################################################################\n");
0459   slhah.addLine("## PARAM_CARD AUTOMATICALY GENERATED BY MG5 FOLLOWING UFO MODEL   ####\n");
0460   slhah.addLine("######################################################################\n");
0461   slhah.addLine("##                                                                  ##\n");
0462   slhah.addLine("##  Width set on Auto will be computed following the information    ##\n");
0463   slhah.addLine("##        present in the decay.py files of the model.               ##\n");
0464   slhah.addLine("##        See  arXiv:1402.1178 for more details.                    ##\n");
0465   slhah.addLine("##                                                                  ##\n");
0466   slhah.addLine("######################################################################\n");
0467   slhah.addLine("\n");
0468   slhah.addLine("###################################\n");
0469   slhah.addLine("## INFORMATION FOR MASS\n");
0470   slhah.addLine("###################################\n");
0471   slhah.addLine("Block mass \n");
0472   slhah.addLine("    6 1.730000e+02 # MT \n");
0473   slhah.addLine("   15 1.777000e+00 # MTA \n");
0474   slhah.addLine("   23 9.118800e+01 # MZ \n");
0475   slhah.addLine("   25 1.250000e+02 # MH \n");
0476   slhah.addLine("## Dependent parameters, given by model restrictions.\n");
0477   slhah.addLine("## Those values should be edited following the \n");
0478   slhah.addLine("## analytical expression. MG5 ignores those values \n");
0479   slhah.addLine("## but they are important for interfacing the output of MG5\n");
0480   slhah.addLine("## to external program such as Pythia.\n");
0481   slhah.addLine("  1 0.000000 # d : 0.0 \n");
0482   slhah.addLine("  2 0.000000 # u : 0.0 \n");
0483   slhah.addLine("  3 0.000000 # s : 0.0 \n");
0484   slhah.addLine("  4 0.000000 # c : 0.0 \n");
0485   slhah.addLine("  5 0.000000 # b : 0.0 \n");
0486   slhah.addLine("  11 0.000000 # e- : 0.0 \n");
0487   slhah.addLine("  12 0.000000 # ve : 0.0 \n");
0488   slhah.addLine("  13 0.000000 # mu- : 0.0 \n");
0489   slhah.addLine("  14 0.000000 # vm : 0.0 \n");
0490   slhah.addLine("  16 0.000000 # vt : 0.0 \n");
0491   slhah.addLine("  21 0.000000 # g : 0.0 \n");
0492   slhah.addLine("  22 0.000000 # a : 0.0 \n");
0493   slhah.addLine(
0494       "  24 80.419002 # w+ : cmath.sqrt(MZ__exp__2/2. + cmath.sqrt(MZ__exp__4/4. - "
0495       "(aEW*cmath.pi*MZ__exp__2)/(Gf*sqrt__2))) \n");
0496   slhah.addLine("\n");
0497   slhah.addLine("###################################\n");
0498   slhah.addLine("## INFORMATION FOR SMINPUTS\n");
0499   slhah.addLine("###################################\n");
0500   slhah.addLine("Block sminputs \n");
0501   slhah.addLine("    1 1.325070e+02 # aEWM1 \n");
0502   slhah.addLine("    2 1.166390e-05 # Gf \n");
0503   slhah.addLine("    3 1.180000e-01 # aS \n");
0504   slhah.addLine("\n");
0505   slhah.addLine("###################################\n");
0506   slhah.addLine("## INFORMATION FOR WOLFENSTEIN\n");
0507   slhah.addLine("###################################\n");
0508   slhah.addLine("Block wolfenstein \n");
0509   slhah.addLine("    1 2.253000e-01 # lamWS \n");
0510   slhah.addLine("    2 8.080000e-01 # AWS \n");
0511   slhah.addLine("    3 1.320000e-01 # rhoWS \n");
0512   slhah.addLine("    4 3.410000e-01 # etaWS \n");
0513   slhah.addLine("\n");
0514   slhah.addLine("###################################\n");
0515   slhah.addLine("## INFORMATION FOR YUKAWA\n");
0516   slhah.addLine("###################################\n");
0517   slhah.addLine("Block yukawa \n");
0518   slhah.addLine("    6 1.730000e+02 # ymt \n");
0519   slhah.addLine("   15 1.777000e+00 # ymtau \n");
0520   slhah.addLine("\n");
0521   slhah.addLine("###################################\n");
0522   slhah.addLine("## INFORMATION FOR DECAY\n");
0523   slhah.addLine("###################################\n");
0524   slhah.addLine("DECAY   6 1.491500e+00 # WT \n");
0525   slhah.addLine("DECAY  15 2.270000e-12 # WTau \n");
0526   slhah.addLine("DECAY  23 2.441404e+00 # WZ \n");
0527   slhah.addLine("DECAY  24 2.047600e+00 # WW \n");
0528   slhah.addLine("DECAY  25 6.382339e-03 # WH \n");
0529   slhah.addLine("## Dependent parameters, given by model restrictions.\n");
0530   slhah.addLine("## Those values should be edited following the \n");
0531   slhah.addLine("## analytical expression. MG5 ignores those values \n");
0532   slhah.addLine("## but they are important for interfacing the output of MG5\n");
0533   slhah.addLine("## to external program such as Pythia.\n");
0534   slhah.addLine("DECAY  1 0.000000 # d : 0.0 \n");
0535   slhah.addLine("DECAY  2 0.000000 # u : 0.0 \n");
0536   slhah.addLine("DECAY  3 0.000000 # s : 0.0 \n");
0537   slhah.addLine("DECAY  4 0.000000 # c : 0.0 \n");
0538   slhah.addLine("DECAY  5 0.000000 # b : 0.0 \n");
0539   slhah.addLine("DECAY  11 0.000000 # e- : 0.0 \n");
0540   slhah.addLine("DECAY  12 0.000000 # ve : 0.0 \n");
0541   slhah.addLine("DECAY  13 0.000000 # mu- : 0.0 \n");
0542   slhah.addLine("DECAY  14 0.000000 # vm : 0.0 \n");
0543   slhah.addLine("DECAY  16 0.000000 # vt : 0.0 \n");
0544   slhah.addLine("DECAY  21 0.000000 # g : 0.0 \n");
0545   slhah.addLine("DECAY  22 0.000000 # a : 0.0\n");
0546 
0547   return slhah;
0548 }
0549 
0550 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0551 void EmbeddingLHEProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0552   //The following says we do not know what parameters are allowed so do no validation
0553   // Please change this to state exactly what you do use, even if it is no parameters
0554   edm::ParameterSetDescription desc;
0555   desc.setUnknown();
0556   descriptions.addDefault(desc);
0557 }
0558 
0559 //define this as a plug-in
0560 DEFINE_FWK_MODULE(EmbeddingLHEProducer);