CMS 3D CMS Logo

GenParticles2HepMCConverter.cc
Go to the documentation of this file.
3 
9 
13 
21 
22 #include <iostream>
23 #include <map>
24 
25 using namespace std;
26 
28 {
29 public:
32 
33  void produce(edm::Event& event, const edm::EventSetup& eventSetup) override;
34 
35 private:
39 
40  std::vector<int> signalParticlePdgIds_;
41 
42 private:
43  inline HepMC::FourVector FourVector(const reco::Candidate::Point& point)
44  {
45  return HepMC::FourVector(10*point.x(), 10*point.y(), 10*point.z(), 0);
46  };
47 
48  inline HepMC::FourVector FourVector(const reco::Candidate::LorentzVector& lvec)
49  {
50  // Avoid negative mass, set minimum m^2 = 0
51  return HepMC::FourVector(lvec.px(), lvec.py(), lvec.pz(), std::hypot(lvec.P(), std::max(0., lvec.mass())));
52  };
53 
54 
55 };
56 
58 {
59  genParticlesToken_ = consumes<reco::CandidateView>(pset.getParameter<edm::InputTag>("genParticles"));
60  genEventInfoToken_ = consumes<GenEventInfoProduct>(pset.getParameter<edm::InputTag>("genEventInfo"));
61  signalParticlePdgIds_ = pset.getParameter<std::vector<int>>("signalParticlePdgIds");
62 
63  produces<edm::HepMCProduct>("unsmeared");
64 }
65 
67 {
68  edm::Handle<reco::CandidateView> genParticlesHandle;
69  event.getByToken(genParticlesToken_, genParticlesHandle);
70 
71  edm::Handle<GenEventInfoProduct> genEventInfoHandle;
72  event.getByToken(genEventInfoToken_, genEventInfoHandle);
73 
74  eventSetup.getData(pTable_);
75 
76  HepMC::GenEvent* hepmc_event = new HepMC::GenEvent();
77  hepmc_event->set_event_number(event.id().event());
78  hepmc_event->set_signal_process_id(genEventInfoHandle->signalProcessID());
79  hepmc_event->set_event_scale(genEventInfoHandle->qScale());
80  hepmc_event->set_alphaQED(genEventInfoHandle->alphaQED());
81  hepmc_event->set_alphaQCD(genEventInfoHandle->alphaQCD());
82 
83  hepmc_event->weights() = genEventInfoHandle->weights();
84 
85  // Set PDF
86  const gen::PdfInfo* pdf = genEventInfoHandle->pdf();
87  const int pdf_id1 = pdf->id.first, pdf_id2 = pdf->id.second;
88  const double pdf_x1 = pdf->x.first, pdf_x2 = pdf->x.second;
89  const double pdf_scalePDF = pdf->scalePDF;
90  const double pdf_xPDF1 = pdf->xPDF.first, pdf_xPDF2 = pdf->xPDF.second;
91  HepMC::PdfInfo hepmc_pdfInfo(pdf_id1, pdf_id2, pdf_x1, pdf_x2, pdf_scalePDF, pdf_xPDF1, pdf_xPDF2);
92  hepmc_event->set_pdf_info(hepmc_pdfInfo);
93 
94  // Prepare list of HepMC::GenParticles
95  std::map<const reco::Candidate*, HepMC::GenParticle*> genCandToHepMCMap;
96  std::vector<HepMC::GenParticle*> hepmc_particles;
97  for ( unsigned int i=0, n=genParticlesHandle->size(); i<n; ++i )
98  {
99  const reco::Candidate* p = &genParticlesHandle->at(i);
100  HepMC::GenParticle* hepmc_particle = new HepMC::GenParticle(FourVector(p->p4()), p->pdgId(), p->status());
101  hepmc_particle->suggest_barcode(i+1);
102 
103  // Assign particle's generated mass from the standard particle data table
104  double particleMass;
105  if ( pTable_->particle(p->pdgId()) ) particleMass = pTable_->particle(p->pdgId())->mass();
106  else particleMass = p->mass();
107 
108  hepmc_particle->set_generated_mass(particleMass);
109 
110  hepmc_particles.push_back(hepmc_particle);
111  genCandToHepMCMap[p] = hepmc_particle;
112  }
113 
114  // Put incident beam particles : proton -> parton vertex
115  const reco::Candidate* parton1 = genParticlesHandle->at(0).daughter(0);
116  const reco::Candidate* parton2 = genParticlesHandle->at(1).daughter(0);
117  HepMC::GenVertex* vertex1 = new HepMC::GenVertex(FourVector(parton1->vertex()));
118  HepMC::GenVertex* vertex2 = new HepMC::GenVertex(FourVector(parton2->vertex()));
119  hepmc_event->add_vertex(vertex1);
120  hepmc_event->add_vertex(vertex2);
121  vertex1->add_particle_in(hepmc_particles[0]);
122  vertex2->add_particle_in(hepmc_particles[1]);
123  hepmc_event->set_beam_particles(hepmc_particles[0], hepmc_particles[1]);
124 
125  // Prepare vertex list
126  typedef std::map<const reco::Candidate*, HepMC::GenVertex*> ParticleToVertexMap;
127  ParticleToVertexMap particleToVertexMap;
128  particleToVertexMap[parton1] = vertex1;
129  particleToVertexMap[parton2] = vertex2;
130  for ( unsigned int i=2, n=genParticlesHandle->size(); i<n; ++i )
131  {
132  const reco::Candidate* p = &genParticlesHandle->at(i);
133 
134  // Connect mother-daughters for the other cases
135  for ( unsigned int j=0, nMothers=p->numberOfMothers(); j<nMothers; ++j )
136  {
137  // Mother-daughter hierarchy defines vertex
138  const reco::Candidate* elder = p->mother(j)->daughter(0);
139  HepMC::GenVertex* vertex;
140  if ( particleToVertexMap.find(elder) == particleToVertexMap.end() )
141  {
142  vertex = new HepMC::GenVertex(FourVector(elder->vertex()));
143  hepmc_event->add_vertex(vertex);
144  particleToVertexMap[elder] = vertex;
145  }
146  else
147  {
148  vertex = particleToVertexMap[elder];
149  }
150 
151  // Vertex is found. Now connect each other
152  const reco::Candidate* mother = p->mother(j);
153  vertex->add_particle_in(genCandToHepMCMap[mother]);
154  vertex->add_particle_out(hepmc_particles[i]);
155  }
156  }
157 
158  // Finalize HepMC event record
159  bool hasSignalVertex = false;
160  if ( !signalParticlePdgIds_.empty() ) {
161  // Loop over all vertices to assign the signal vertex, decaying to a signal particle
162  for ( auto v = hepmc_event->vertices_begin(); v != hepmc_event->vertices_end(); ++v ) {
163  for ( auto p = (*v)->particles_begin(HepMC::children);
164  p != (*v)->particles_end(HepMC::children); ++p ) {
165  const int pdgId = (*p)->pdg_id();
166  if ( std::find(signalParticlePdgIds_.begin(), signalParticlePdgIds_.end(), pdgId) != signalParticlePdgIds_.end() ) {
167  hepmc_event->set_signal_process_vertex(*v);
168  hasSignalVertex = true;
169  break;
170  }
171  }
172  if ( hasSignalVertex ) break;
173  }
174  }
175  // Set the default signal vertex if still not set
176  if ( !hasSignalVertex ) hepmc_event->set_signal_process_vertex(*(vertex1->vertices_begin()));
177 
178  std::unique_ptr<edm::HepMCProduct> hepmc_product(new edm::HepMCProduct());
179  hepmc_product->addHepMCData(hepmc_event);
180  event.put(std::move(hepmc_product), "unsmeared");
181 
182 }
183 
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoToken_
const PDF * pdf() const
double alphaQCD() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
std::pair< double, double > x
Definition: PdfInfo.h:13
edm::ESHandle< ParticleDataTable > pTable_
size_type size() const
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
HepMC::FourVector FourVector(const reco::Candidate::Point &point)
edm::EDGetTokenT< reco::CandidateView > genParticlesToken_
void getData(T &iHolder) const
Definition: EventSetup.h:78
std::pair< double, double > xPDF
Definition: PdfInfo.h:14
virtual int status() const =0
status word
double qScale() const
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
HepMC::FourVector FourVector(const reco::Candidate::LorentzVector &lvec)
virtual int pdgId() const =0
PDG identifier.
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
double alphaQED() const
std::pair< int, int > id
Definition: PdfInfo.h:12
virtual double mass() const =0
mass
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
edm::EventID id() const
Definition: EventBase.h:60
std::vector< double > & weights()
const_reference at(size_type pos) const
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
virtual const Point & vertex() const =0
vertex position
unsigned int signalProcessID() const
GenParticles2HepMCConverter(const edm::ParameterSet &pset)
double scalePDF
Definition: PdfInfo.h:15
def move(src, dest)
Definition: eostools.py:510
*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: event.py:1