CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/GeneratorInterface/RivetInterface/plugins/GenParticles2HepMCConverter.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/Frameworkfwd.h"
00002 #include "FWCore/Framework/interface/EDProducer.h"
00003 
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/Framework/interface/EventSetup.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/Run.h"
00008 #include "FWCore/Framework/interface/MakerMacros.h"
00009 
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012 #include "DataFormats/Common/interface/Handle.h"
00013 
00014 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00015 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00016 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00017 //#include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
00018 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00019 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00020 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00021 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00022 
00023 #include <iostream>
00024 #include <map>
00025 
00026 using namespace std;
00027 
00028 class GenParticles2HepMCConverter : public edm::EDProducer
00029 {
00030 public:
00031   explicit GenParticles2HepMCConverter(const edm::ParameterSet& pset);
00032   ~GenParticles2HepMCConverter() {};
00033 
00034   void beginRun(edm::Run& run, const edm::EventSetup& eventSetup);
00035   void produce(edm::Event& event, const edm::EventSetup& eventSetup);
00036 
00037 private:
00038 //  edm::InputTag lheEventLabel_;
00039   edm::InputTag genParticlesLabel_;
00040 //  edm::InputTag genRunInfoLabel_;
00041   edm::InputTag genEventInfoLabel_;
00042   edm::ESHandle<ParticleDataTable> pTable_;
00043 
00044 private:
00045   inline HepMC::FourVector FourVector(const reco::Candidate::Point& point)
00046   {
00047     return HepMC::FourVector(10*point.x(), 10*point.y(), 10*point.z(), 0);
00048   };
00049 
00050   inline HepMC::FourVector FourVector(const reco::Candidate::LorentzVector& lvec)
00051   {
00052     // Avoid negative mass, set minimum m^2 = 0
00053     return HepMC::FourVector(lvec.px(), lvec.py(), lvec.pz(), std::hypot(lvec.P(), std::max(0., lvec.mass())));
00054   };
00055 
00056 
00057 };
00058 
00059 GenParticles2HepMCConverter::GenParticles2HepMCConverter(const edm::ParameterSet& pset)
00060 {
00061 //  lheEventLabel_ = pset.getParameter<edm::InputTag>("lheEvent");
00062   genParticlesLabel_ = pset.getParameter<edm::InputTag>("genParticles");
00063   //genRunInfoLabel_ = pset.getParameter<edm::InputTag>("genRunInfo");
00064   genEventInfoLabel_ = pset.getParameter<edm::InputTag>("genEventInfo");
00065 
00066   produces<edm::HepMCProduct>();
00067 }
00068 
00069 void GenParticles2HepMCConverter::beginRun(edm::Run& run, const edm::EventSetup& eventSetup)
00070 {
00071   //edm::Handle<GenRunInfoProduct> genRunInfoHandle;
00072   //event.getByLabel(genRunInfoLabel_, genRunInfoHandle);
00073   // const double xsecIn = genRunInfoHandle->internalXSec().value();
00074   // const double xsecInErr = genRunInfoHandle->internalXSec().error();
00075   // const double xsecLO = genRunInfoHandle->externalXSecLO().value();
00076   // const double xsecLOErr = genRunInfoHandle->externalXSecLO().error();
00077   // const double xsecNLO = genRunInfoHandle->externalXSecNLO().value();
00078   // const double xsecNLOErr = genRunInfoHandle->externalXSecNLO().error();
00079 }
00080 
00081 void GenParticles2HepMCConverter::produce(edm::Event& event, const edm::EventSetup& eventSetup)
00082 {
00083 //  edm::Handle<LHEEventProduct> lheEventHandle;
00084 //  event.getByLabel(lheEventLabel_, lheEventHandle);
00085 
00086   edm::Handle<reco::CandidateView> genParticlesHandle;
00087   event.getByLabel(genParticlesLabel_, genParticlesHandle);
00088 
00089   edm::Handle<GenEventInfoProduct> genEventInfoHandle;
00090   event.getByLabel(genEventInfoLabel_, genEventInfoHandle);
00091 
00092   eventSetup.getData(pTable_);
00093 
00094   HepMC::GenEvent* hepmc_event = new HepMC::GenEvent();
00095   hepmc_event->set_event_number(event.id().event());
00096   hepmc_event->set_signal_process_id(genEventInfoHandle->signalProcessID());
00097   hepmc_event->set_event_scale(genEventInfoHandle->qScale());
00098   hepmc_event->set_alphaQED(genEventInfoHandle->alphaQED());
00099   hepmc_event->set_alphaQCD(genEventInfoHandle->alphaQCD());
00100 
00101   hepmc_event->weights() = genEventInfoHandle->weights();
00102 
00103   // Set PDF
00104   const gen::PdfInfo* pdf = genEventInfoHandle->pdf();
00105   const int pdf_id1 = pdf->id.first, pdf_id2 = pdf->id.second;
00106   const double pdf_x1 = pdf->x.first, pdf_x2 = pdf->x.second;
00107   const double pdf_scalePDF = pdf->scalePDF;
00108   const double pdf_xPDF1 = pdf->xPDF.first, pdf_xPDF2 = pdf->xPDF.second;
00109   HepMC::PdfInfo hepmc_pdfInfo(pdf_id1, pdf_id2, pdf_x1, pdf_x2, pdf_scalePDF, pdf_xPDF1, pdf_xPDF2);
00110   hepmc_event->set_pdf_info(hepmc_pdfInfo);
00111 
00112   // Load LHE
00113 //  const lhef::HEPEUP& lheEvent = lheEventHandle->hepeup();
00114 //  std::vector<int> lhe_meIndex; // Particle indices with preserved mass, status=2
00115 //  for ( int i=0, n=lheEvent.ISTUP.size(); i<n; ++i )
00116 //  {
00117 //    if ( lheEvent.ISTUP[i] == 2 ) lhe_meIndex.push_back(i);
00118 //  }
00119 
00120   // Prepare list of HepMC::GenParticles
00121   std::map<const reco::Candidate*, HepMC::GenParticle*> genCandToHepMCMap;
00122   std::vector<HepMC::GenParticle*> hepmc_particles;
00123   for ( unsigned int i=0, n=genParticlesHandle->size(); i<n; ++i )
00124   {
00125     const reco::Candidate* p = &genParticlesHandle->at(i);
00126     HepMC::GenParticle* hepmc_particle = new HepMC::GenParticle(FourVector(p->p4()), p->pdgId(), p->status());
00127     hepmc_particle->suggest_barcode(i+1);
00128 
00129     // Assign particle's generated mass from the standard particle data table
00130     double particleMass = pTable_->particle(p->pdgId())->mass();
00131 //    // Re-assign generated mass from LHE, find particle among the LHE
00132 //    for ( unsigned int j=0, m=lhe_meIndex.size(); j<m; ++j )
00133 //    {
00134 //      const unsigned int lheIndex = lhe_meIndex[j];
00135 //      if ( p->pdgId() != lheEvent.IDUP[lheIndex] ) continue;
00136 //
00137 //      const lhef::HEPEUP::FiveVector& vp = lheEvent.PUP[lheIndex];
00138 //      if ( std::abs(vp[0] - p->px()) > 1e-7 or std::abs(vp[1] - p->py()) > 1e-7 ) continue;
00139 //      if ( std::abs(vp[2] - p->pz()) > 1e-7 or std::abs(vp[3] - p->energy()) > 1e-7 ) continue;
00140 //
00141 //      particleMass = vp[4];
00142 //      break;
00143 //    }
00144     hepmc_particle->set_generated_mass(particleMass);
00145 
00146     hepmc_particles.push_back(hepmc_particle);
00147     genCandToHepMCMap[p] = hepmc_particle;
00148   }
00149 
00150   // Put incident beam particles : proton -> parton vertex
00151   const reco::Candidate* parton1 = genParticlesHandle->at(0).daughter(0);
00152   const reco::Candidate* parton2 = genParticlesHandle->at(1).daughter(0);
00153   HepMC::GenVertex* vertex1 = new HepMC::GenVertex(FourVector(parton1->vertex()));
00154   HepMC::GenVertex* vertex2 = new HepMC::GenVertex(FourVector(parton2->vertex()));
00155   hepmc_event->add_vertex(vertex1);
00156   hepmc_event->add_vertex(vertex2);
00157   //hepmc_particles[0]->set_status(4);
00158   //hepmc_particles[1]->set_status(4);
00159   vertex1->add_particle_in(hepmc_particles[0]);
00160   vertex2->add_particle_in(hepmc_particles[1]);
00161   hepmc_event->set_beam_particles(hepmc_particles[0], hepmc_particles[1]);
00162 
00163   // Prepare vertex list
00164   typedef std::map<const reco::Candidate*, HepMC::GenVertex*> ParticleToVertexMap;
00165   ParticleToVertexMap particleToVertexMap;
00166   particleToVertexMap[parton1] = vertex1;
00167   particleToVertexMap[parton2] = vertex2;
00168   for ( unsigned int i=2, n=genParticlesHandle->size(); i<n; ++i )
00169   {
00170     const reco::Candidate* p = &genParticlesHandle->at(i);
00171 
00172     // Connect mother-daughters for the other cases
00173     for ( unsigned int j=0, nMothers=p->numberOfMothers(); j<nMothers; ++j )
00174     {
00175       // Mother-daughter hierarchy defines vertex
00176       const reco::Candidate* elder = p->mother(j)->daughter(0);
00177       HepMC::GenVertex* vertex;
00178       if ( particleToVertexMap.find(elder) == particleToVertexMap.end() )
00179       {
00180         vertex = new HepMC::GenVertex(FourVector(elder->vertex()));
00181         hepmc_event->add_vertex(vertex);
00182         particleToVertexMap[elder] = vertex;
00183       }
00184       else
00185       {
00186         vertex = particleToVertexMap[elder];
00187       }
00188 
00189       // Vertex is found. Now connect each other
00190       const reco::Candidate* mother = p->mother(j);
00191       vertex->add_particle_in(genCandToHepMCMap[mother]);
00192       vertex->add_particle_out(hepmc_particles[i]);
00193     }
00194   }
00195 
00196   // Finalize HepMC event record
00197   hepmc_event->set_signal_process_vertex(*(vertex1->vertices_begin()));
00198 
00199   std::auto_ptr<edm::HepMCProduct> hepmc_product(new edm::HepMCProduct());
00200   hepmc_product->addHepMCData(hepmc_event);
00201   event.put(hepmc_product);
00202 
00203 }
00204 
00205 DEFINE_FWK_MODULE(GenParticles2HepMCConverter);