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  HepMC::GenParticle* hepmc_parton1 = nullptr, * hepmc_parton2 = nullptr;
97  std::vector<HepMC::GenParticle*> hepmc_particles;
98  const reco::Candidate* parton1 = nullptr, * parton2 = nullptr;
99  for ( unsigned int i=0, n=genParticlesHandle->size(); i<n; ++i )
100  {
101  const reco::Candidate* p = &genParticlesHandle->at(i);
102  HepMC::GenParticle* hepmc_particle = new HepMC::GenParticle(FourVector(p->p4()), p->pdgId(), p->status());
103  hepmc_particle->suggest_barcode(i+1);
104 
105  // Assign particle's generated mass from the standard particle data table
106  double particleMass;
107  if ( pTable_->particle(p->pdgId()) ) particleMass = pTable_->particle(p->pdgId())->mass();
108  else particleMass = p->mass();
109 
110  hepmc_particle->set_generated_mass(particleMass);
111 
112  hepmc_particles.push_back(hepmc_particle);
113  genCandToHepMCMap[p] = hepmc_particle;
114 
115  // Find incident proton pair
116  if ( p->mother() == nullptr and std::abs(p->eta()) > 5 and std::abs(p->pz()) > 1000 ) {
117  if ( !parton1 and p->pz() > 0 ) {
118  parton1 = p;
119  hepmc_parton1 = hepmc_particle;
120  }
121  else if ( !parton2 and p->pz() < 0 ) {
122  parton2 = p;
123  hepmc_parton2 = hepmc_particle;
124  }
125  }
126  }
127 
128  // Put incident beam particles : proton -> parton vertex
129  HepMC::GenVertex* vertex1 = new HepMC::GenVertex(FourVector(parton1->vertex()));
130  HepMC::GenVertex* vertex2 = new HepMC::GenVertex(FourVector(parton2->vertex()));
131  hepmc_event->add_vertex(vertex1);
132  hepmc_event->add_vertex(vertex2);
133  vertex1->add_particle_in(hepmc_parton1);
134  vertex2->add_particle_in(hepmc_parton2);
135  //hepmc_event->set_beam_particles(hepmc_parton1, hepmc_parton2);
136 
137  // Prepare vertex list
138  typedef std::map<const reco::Candidate*, HepMC::GenVertex*> ParticleToVertexMap;
139  ParticleToVertexMap particleToVertexMap;
140  particleToVertexMap[parton1] = vertex1;
141  particleToVertexMap[parton2] = vertex2;
142  for ( unsigned int i=0, n=genParticlesHandle->size(); i<n; ++i )
143  {
144  const reco::Candidate* p = &genParticlesHandle->at(i);
145  if ( p == parton1 or p == parton2 ) continue;
146 
147  // Connect mother-daughters for the other cases
148  for ( unsigned int j=0, nMothers=p->numberOfMothers(); j<nMothers; ++j )
149  {
150  // Mother-daughter hierarchy defines vertex
151  const reco::Candidate* elder = p->mother(j)->daughter(0);
152  HepMC::GenVertex* vertex;
153  if ( particleToVertexMap.find(elder) == particleToVertexMap.end() )
154  {
155  vertex = new HepMC::GenVertex(FourVector(elder->vertex()));
156  hepmc_event->add_vertex(vertex);
157  particleToVertexMap[elder] = vertex;
158  }
159  else
160  {
161  vertex = particleToVertexMap[elder];
162  }
163 
164  // Vertex is found. Now connect each other
165  const reco::Candidate* mother = p->mother(j);
166  vertex->add_particle_in(genCandToHepMCMap[mother]);
167  vertex->add_particle_out(hepmc_particles[i]);
168  }
169  }
170 
171  // Finalize HepMC event record
172  bool hasSignalVertex = false;
173  if ( !signalParticlePdgIds_.empty() ) {
174  // Loop over all vertices to assign the signal vertex, decaying to a signal particle
175  for ( auto v = hepmc_event->vertices_begin(); v != hepmc_event->vertices_end(); ++v ) {
176  for ( auto p = (*v)->particles_begin(HepMC::children);
177  p != (*v)->particles_end(HepMC::children); ++p ) {
178  const int pdgId = (*p)->pdg_id();
179  if ( std::find(signalParticlePdgIds_.begin(), signalParticlePdgIds_.end(), pdgId) != signalParticlePdgIds_.end() ) {
180  hepmc_event->set_signal_process_vertex(*v);
181  hasSignalVertex = true;
182  break;
183  }
184  }
185  if ( hasSignalVertex ) break;
186  }
187  }
188  // Set the default signal vertex if still not set
189  if ( !hasSignalVertex ) hepmc_event->set_signal_process_vertex(*(vertex1->vertices_begin()));
190 
191  std::unique_ptr<edm::HepMCProduct> hepmc_product(new edm::HepMCProduct());
192  hepmc_product->addHepMCData(hepmc_event);
193  event.put(std::move(hepmc_product), "unsmeared");
194 
195 }
196 
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
virtual double pz() const =0
z coordinate of momentum vector
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:81
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.
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double alphaQED() const
std::pair< int, int > id
Definition: PdfInfo.h:12
virtual double eta() const =0
momentum pseudorapidity
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