CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenParticles2HepMCConverter.cc
Go to the documentation of this file.
3 
9 
13 
17 //#include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
22 
23 #include <iostream>
24 #include <map>
25 
26 using namespace std;
27 
29 {
30 public:
31  explicit GenParticles2HepMCConverter(const edm::ParameterSet& pset);
33 
34  void beginRun(edm::Run& run, const edm::EventSetup& eventSetup);
35  void produce(edm::Event& event, const edm::EventSetup& eventSetup);
36 
37 private:
38 // edm::InputTag lheEventLabel_;
40 // edm::InputTag genRunInfoLabel_;
43 
44 private:
45  inline HepMC::FourVector FourVector(const reco::Candidate::Point& point)
46  {
47  return HepMC::FourVector(10*point.x(), 10*point.y(), 10*point.z(), 0);
48  };
49 
50  inline HepMC::FourVector FourVector(const reco::Candidate::LorentzVector& lvec)
51  {
52  // Avoid negative mass, set minimum m^2 = 0
53  return HepMC::FourVector(lvec.px(), lvec.py(), lvec.pz(), std::hypot(lvec.P(), std::max(0., lvec.mass())));
54  };
55 
56 
57 };
58 
60 {
61 // lheEventLabel_ = pset.getParameter<edm::InputTag>("lheEvent");
62  genParticlesLabel_ = pset.getParameter<edm::InputTag>("genParticles");
63  //genRunInfoLabel_ = pset.getParameter<edm::InputTag>("genRunInfo");
64  genEventInfoLabel_ = pset.getParameter<edm::InputTag>("genEventInfo");
65 
66  produces<edm::HepMCProduct>();
67 }
68 
70 {
71  //edm::Handle<GenRunInfoProduct> genRunInfoHandle;
72  //event.getByLabel(genRunInfoLabel_, genRunInfoHandle);
73  // const double xsecIn = genRunInfoHandle->internalXSec().value();
74  // const double xsecInErr = genRunInfoHandle->internalXSec().error();
75  // const double xsecLO = genRunInfoHandle->externalXSecLO().value();
76  // const double xsecLOErr = genRunInfoHandle->externalXSecLO().error();
77  // const double xsecNLO = genRunInfoHandle->externalXSecNLO().value();
78  // const double xsecNLOErr = genRunInfoHandle->externalXSecNLO().error();
79 }
80 
82 {
83 // edm::Handle<LHEEventProduct> lheEventHandle;
84 // event.getByLabel(lheEventLabel_, lheEventHandle);
85 
86  edm::Handle<reco::CandidateView> genParticlesHandle;
87  event.getByLabel(genParticlesLabel_, genParticlesHandle);
88 
89  edm::Handle<GenEventInfoProduct> genEventInfoHandle;
90  event.getByLabel(genEventInfoLabel_, genEventInfoHandle);
91 
92  eventSetup.getData(pTable_);
93 
94  HepMC::GenEvent* hepmc_event = new HepMC::GenEvent();
95  hepmc_event->set_event_number(event.id().event());
96  hepmc_event->set_signal_process_id(genEventInfoHandle->signalProcessID());
97  hepmc_event->set_event_scale(genEventInfoHandle->qScale());
98  hepmc_event->set_alphaQED(genEventInfoHandle->alphaQED());
99  hepmc_event->set_alphaQCD(genEventInfoHandle->alphaQCD());
100 
101  hepmc_event->weights() = genEventInfoHandle->weights();
102 
103  // Set PDF
104  const gen::PdfInfo* pdf = genEventInfoHandle->pdf();
105  const int pdf_id1 = pdf->id.first, pdf_id2 = pdf->id.second;
106  const double pdf_x1 = pdf->x.first, pdf_x2 = pdf->x.second;
107  const double pdf_scalePDF = pdf->scalePDF;
108  const double pdf_xPDF1 = pdf->xPDF.first, pdf_xPDF2 = pdf->xPDF.second;
109  HepMC::PdfInfo hepmc_pdfInfo(pdf_id1, pdf_id2, pdf_x1, pdf_x2, pdf_scalePDF, pdf_xPDF1, pdf_xPDF2);
110  hepmc_event->set_pdf_info(hepmc_pdfInfo);
111 
112  // Load LHE
113 // const lhef::HEPEUP& lheEvent = lheEventHandle->hepeup();
114 // std::vector<int> lhe_meIndex; // Particle indices with preserved mass, status=2
115 // for ( int i=0, n=lheEvent.ISTUP.size(); i<n; ++i )
116 // {
117 // if ( lheEvent.ISTUP[i] == 2 ) lhe_meIndex.push_back(i);
118 // }
119 
120  // Prepare list of HepMC::GenParticles
121  std::map<const reco::Candidate*, HepMC::GenParticle*> genCandToHepMCMap;
122  std::vector<HepMC::GenParticle*> hepmc_particles;
123  for ( unsigned int i=0, n=genParticlesHandle->size(); i<n; ++i )
124  {
125  const reco::Candidate* p = &genParticlesHandle->at(i);
126  HepMC::GenParticle* hepmc_particle = new HepMC::GenParticle(FourVector(p->p4()), p->pdgId(), p->status());
127  hepmc_particle->suggest_barcode(i+1);
128 
129  // Assign particle's generated mass from the standard particle data table
130  double particleMass;
131  if ( pTable_->particle(p->pdgId()) ) particleMass = pTable_->particle(p->pdgId())->mass();
132  else particleMass = p->mass();
133 // // Re-assign generated mass from LHE, find particle among the LHE
134 // for ( unsigned int j=0, m=lhe_meIndex.size(); j<m; ++j )
135 // {
136 // const unsigned int lheIndex = lhe_meIndex[j];
137 // if ( p->pdgId() != lheEvent.IDUP[lheIndex] ) continue;
138 //
139 // const lhef::HEPEUP::FiveVector& vp = lheEvent.PUP[lheIndex];
140 // if ( std::abs(vp[0] - p->px()) > 1e-7 or std::abs(vp[1] - p->py()) > 1e-7 ) continue;
141 // if ( std::abs(vp[2] - p->pz()) > 1e-7 or std::abs(vp[3] - p->energy()) > 1e-7 ) continue;
142 //
143 // particleMass = vp[4];
144 // break;
145 // }
146  hepmc_particle->set_generated_mass(particleMass);
147 
148  hepmc_particles.push_back(hepmc_particle);
149  genCandToHepMCMap[p] = hepmc_particle;
150  }
151 
152  // Put incident beam particles : proton -> parton vertex
153  const reco::Candidate* parton1 = genParticlesHandle->at(0).daughter(0);
154  const reco::Candidate* parton2 = genParticlesHandle->at(1).daughter(0);
155  HepMC::GenVertex* vertex1 = new HepMC::GenVertex(FourVector(parton1->vertex()));
156  HepMC::GenVertex* vertex2 = new HepMC::GenVertex(FourVector(parton2->vertex()));
157  hepmc_event->add_vertex(vertex1);
158  hepmc_event->add_vertex(vertex2);
159  //hepmc_particles[0]->set_status(4);
160  //hepmc_particles[1]->set_status(4);
161  vertex1->add_particle_in(hepmc_particles[0]);
162  vertex2->add_particle_in(hepmc_particles[1]);
163  hepmc_event->set_beam_particles(hepmc_particles[0], hepmc_particles[1]);
164 
165  // Prepare vertex list
166  typedef std::map<const reco::Candidate*, HepMC::GenVertex*> ParticleToVertexMap;
167  ParticleToVertexMap particleToVertexMap;
168  particleToVertexMap[parton1] = vertex1;
169  particleToVertexMap[parton2] = vertex2;
170  for ( unsigned int i=2, n=genParticlesHandle->size(); i<n; ++i )
171  {
172  const reco::Candidate* p = &genParticlesHandle->at(i);
173 
174  // Connect mother-daughters for the other cases
175  for ( unsigned int j=0, nMothers=p->numberOfMothers(); j<nMothers; ++j )
176  {
177  // Mother-daughter hierarchy defines vertex
178  const reco::Candidate* elder = p->mother(j)->daughter(0);
179  HepMC::GenVertex* vertex;
180  if ( particleToVertexMap.find(elder) == particleToVertexMap.end() )
181  {
182  vertex = new HepMC::GenVertex(FourVector(elder->vertex()));
183  hepmc_event->add_vertex(vertex);
184  particleToVertexMap[elder] = vertex;
185  }
186  else
187  {
188  vertex = particleToVertexMap[elder];
189  }
190 
191  // Vertex is found. Now connect each other
192  const reco::Candidate* mother = p->mother(j);
193  vertex->add_particle_in(genCandToHepMCMap[mother]);
194  vertex->add_particle_out(hepmc_particles[i]);
195  }
196  }
197 
198  // Finalize HepMC event record
199  hepmc_event->set_signal_process_vertex(*(vertex1->vertices_begin()));
200 
201  std::auto_ptr<edm::HepMCProduct> hepmc_product(new edm::HepMCProduct());
202  hepmc_product->addHepMCData(hepmc_event);
203  event.put(hepmc_product);
204 
205 }
206 
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual float mass() const =0
mass
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
void produce(edm::Event &event, const edm::EventSetup &eventSetup)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual int status() const =0
status word
std::pair< double, double > x
Definition: PdfInfo.h:11
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
edm::ESHandle< ParticleDataTable > pTable_
HepMC::FourVector FourVector(const reco::Candidate::Point &point)
void getData(T &iHolder) const
Definition: EventSetup.h:78
std::pair< double, double > xPDF
Definition: PdfInfo.h:12
const T & max(const T &a, const T &b)
HepMC::FourVector FourVector(const reco::Candidate::LorentzVector &lvec)
void beginRun(edm::Run &run, const edm::EventSetup &eventSetup)
int j
Definition: DBlmapReader.cc:9
virtual const Point & vertex() const =0
vertex position
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
virtual int pdgId() const =0
PDG identifier.
std::pair< int, int > id
Definition: PdfInfo.h:10
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
edm::EventID id() const
Definition: EventBase.h:56
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
GenParticles2HepMCConverter(const edm::ParameterSet &pset)
double scalePDF
Definition: PdfInfo.h:13
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
Definition: Run.h:41
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector