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 public:
31 
32  void beginRun(edm::Run const& iRun, edm::EventSetup const&) override;
33  void produce(edm::Event& event, const edm::EventSetup& eventSetup) override;
34 
35 private:
40 
41  std::vector<int> signalParticlePdgIds_;
42  const double cmEnergy_;
43  HepMC::GenCrossSection xsec_;
44 
45 private:
46  inline HepMC::FourVector FourVector(const reco::Candidate::Point& point) {
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  // Avoid negative mass, set minimum m^2 = 0
52  return HepMC::FourVector(lvec.px(), lvec.py(), lvec.pz(), std::hypot(lvec.P(), std::max(0., lvec.mass())));
53  };
54 };
55 
57  // dummy value to set incident proton pz for particle gun samples
58  : cmEnergy_(pset.getUntrackedParameter<double>("cmEnergy", 13000)) {
59  genParticlesToken_ = consumes<reco::CandidateView>(pset.getParameter<edm::InputTag>("genParticles"));
60  genEventInfoToken_ = consumes<GenEventInfoProduct>(pset.getParameter<edm::InputTag>("genEventInfo"));
61  genRunInfoToken_ = consumes<GenRunInfoProduct, edm::InRun>(pset.getParameter<edm::InputTag>("genEventInfo"));
62  pTable_ = esConsumes<HepPDT::ParticleDataTable, PDTRecord>();
63  signalParticlePdgIds_ = pset.getParameter<std::vector<int>>("signalParticlePdgIds");
64 
65  produces<edm::HepMCProduct>("unsmeared");
66 }
67 
69  edm::Handle<GenRunInfoProduct> genRunInfoHandle;
70  iRun.getByToken(genRunInfoToken_, genRunInfoHandle);
71 
72  xsec_.set_cross_section(genRunInfoHandle->internalXSec().value(), genRunInfoHandle->internalXSec().error());
73 }
74 
76  edm::Handle<reco::CandidateView> genParticlesHandle;
77  event.getByToken(genParticlesToken_, genParticlesHandle);
78 
79  edm::Handle<GenEventInfoProduct> genEventInfoHandle;
80  event.getByToken(genEventInfoToken_, genEventInfoHandle);
81 
82  auto const& pTableData = eventSetup.getData(pTable_);
83 
84  HepMC::GenEvent* hepmc_event = new HepMC::GenEvent();
85  hepmc_event->set_event_number(event.id().event());
86  hepmc_event->set_signal_process_id(genEventInfoHandle->signalProcessID());
87  hepmc_event->set_event_scale(genEventInfoHandle->qScale());
88  hepmc_event->set_alphaQED(genEventInfoHandle->alphaQED());
89  hepmc_event->set_alphaQCD(genEventInfoHandle->alphaQCD());
90 
91  hepmc_event->weights() = genEventInfoHandle->weights();
92 
93  hepmc_event->set_cross_section(xsec_);
94 
95  // Set PDF
96  const gen::PdfInfo* pdf = genEventInfoHandle->pdf();
97  if (pdf != nullptr) {
98  const int pdf_id1 = pdf->id.first, pdf_id2 = pdf->id.second;
99  const double pdf_x1 = pdf->x.first, pdf_x2 = pdf->x.second;
100  const double pdf_scalePDF = pdf->scalePDF;
101  const double pdf_xPDF1 = pdf->xPDF.first, pdf_xPDF2 = pdf->xPDF.second;
102  HepMC::PdfInfo hepmc_pdfInfo(pdf_id1, pdf_id2, pdf_x1, pdf_x2, pdf_scalePDF, pdf_xPDF1, pdf_xPDF2);
103  hepmc_event->set_pdf_info(hepmc_pdfInfo);
104  }
105 
106  // Prepare list of HepMC::GenParticles
107  std::map<const reco::Candidate*, HepMC::GenParticle*> genCandToHepMCMap;
108  HepMC::GenParticle *hepmc_parton1 = nullptr, *hepmc_parton2 = nullptr;
109  std::vector<HepMC::GenParticle*> hepmc_particles;
110  const reco::Candidate *parton1 = nullptr, *parton2 = nullptr;
111  for (unsigned int i = 0, n = genParticlesHandle->size(); i < n; ++i) {
112  const reco::Candidate* p = &genParticlesHandle->at(i);
113  HepMC::GenParticle* hepmc_particle = new HepMC::GenParticle(FourVector(p->p4()), p->pdgId(), p->status());
114  hepmc_particle->suggest_barcode(i + 1);
115 
116  // Assign particle's generated mass from the standard particle data table
117  double particleMass;
118  if (pTableData.particle(p->pdgId()))
119  particleMass = pTableData.particle(p->pdgId())->mass();
120  else
121  particleMass = p->mass();
122 
123  hepmc_particle->set_generated_mass(particleMass);
124 
125  hepmc_particles.push_back(hepmc_particle);
126  genCandToHepMCMap[p] = hepmc_particle;
127 
128  // Find incident proton pair
129  if (p->mother() == nullptr and std::abs(p->eta()) > 5 and std::abs(p->pz()) > 1000) {
130  if (!parton1 and p->pz() > 0) {
131  parton1 = p;
132  hepmc_parton1 = hepmc_particle;
133  } else if (!parton2 and p->pz() < 0) {
134  parton2 = p;
135  hepmc_parton2 = hepmc_particle;
136  }
137  }
138  }
139 
140  HepMC::GenVertex* vertex1 = nullptr;
141  HepMC::GenVertex* vertex2 = nullptr;
142  if (parton1 == nullptr || parton2 == nullptr) {
143  // Particle gun samples do not have incident partons. Put dummy incident particle and prod vertex
144  // Note: leave parton1 and parton2 as nullptr since it is not used anymore after creating hepmc_parton1 and 2
145  const reco::Candidate::LorentzVector nullP4(0, 0, 0, 0);
146  const reco::Candidate::LorentzVector beamP4(0, 0, cmEnergy_ / 2, cmEnergy_ / 2);
147  vertex1 = new HepMC::GenVertex(FourVector(nullP4));
148  vertex2 = new HepMC::GenVertex(FourVector(nullP4));
149  hepmc_parton1 = new HepMC::GenParticle(FourVector(+beamP4), 2212, 4);
150  hepmc_parton2 = new HepMC::GenParticle(FourVector(-beamP4), 2212, 4);
151  } else {
152  // Put incident beam particles : proton -> parton vertex
153  vertex1 = new HepMC::GenVertex(FourVector(parton1->vertex()));
154  vertex2 = new HepMC::GenVertex(FourVector(parton2->vertex()));
155  }
156  hepmc_event->add_vertex(vertex1);
157  hepmc_event->add_vertex(vertex2);
158  vertex1->add_particle_in(hepmc_parton1);
159  vertex2->add_particle_in(hepmc_parton2);
160  //hepmc_event->set_beam_particles(hepmc_parton1, hepmc_parton2);
161 
162  // Prepare vertex list
163  typedef std::map<const reco::Candidate*, HepMC::GenVertex*> ParticleToVertexMap;
164  ParticleToVertexMap particleToVertexMap;
165  particleToVertexMap[parton1] = vertex1;
166  particleToVertexMap[parton2] = vertex2;
167  for (unsigned int i = 0, n = genParticlesHandle->size(); i < n; ++i) {
168  const reco::Candidate* p = &genParticlesHandle->at(i);
169  if (p == parton1 or p == parton2)
170  continue;
171 
172  // Connect mother-daughters for the other cases
173  for (unsigned int j = 0, nMothers = p->numberOfMothers(); j < nMothers; ++j) {
174  // Mother-daughter hierarchy defines vertex
175  const reco::Candidate* elder = p->mother(j)->daughter(0);
176  HepMC::GenVertex* vertex;
177  if (particleToVertexMap.find(elder) == particleToVertexMap.end()) {
178  vertex = new HepMC::GenVertex(FourVector(elder->vertex()));
179  hepmc_event->add_vertex(vertex);
180  particleToVertexMap[elder] = vertex;
181  } else {
182  vertex = particleToVertexMap[elder];
183  }
184 
185  // Vertex is found. Now connect each other
186  const reco::Candidate* mother = p->mother(j);
187  vertex->add_particle_in(genCandToHepMCMap[mother]);
188  vertex->add_particle_out(hepmc_particles[i]);
189  }
190  }
191 
192  // Finalize HepMC event record
193  bool hasSignalVertex = false;
194  if (!signalParticlePdgIds_.empty()) {
195  // Loop over all vertices to assign the signal vertex, decaying to a signal particle
196  for (auto v = hepmc_event->vertices_begin(); v != hepmc_event->vertices_end(); ++v) {
197  for (auto p = (*v)->particles_begin(HepMC::children); p != (*v)->particles_end(HepMC::children); ++p) {
198  const int pdgId = (*p)->pdg_id();
200  signalParticlePdgIds_.end()) {
201  hepmc_event->set_signal_process_vertex(*v);
202  hasSignalVertex = true;
203  break;
204  }
205  }
206  if (hasSignalVertex)
207  break;
208  }
209  }
210  // Set the default signal vertex if still not set
211  if (!hasSignalVertex)
212  hepmc_event->set_signal_process_vertex(*(vertex1->vertices_begin()));
213 
214  std::unique_ptr<edm::HepMCProduct> hepmc_product(new edm::HepMCProduct());
215  hepmc_product->addHepMCData(hepmc_event);
216  event.put(std::move(hepmc_product), "unsmeared");
217 }
218 
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
edm::EDGetTokenT< GenRunInfoProduct > genRunInfoToken_
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoToken_
const XSec & internalXSec() const
void beginRun(edm::Run const &iRun, edm::EventSetup const &) override
std::pair< double, double > x
Definition: PdfInfo.h:13
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
HepMC::FourVector FourVector(const reco::Candidate::Point &point)
edm::EDGetTokenT< reco::CandidateView > genParticlesToken_
size_type size() const
edm::ESGetToken< HepPDT::ParticleDataTable, PDTRecord > pTable_
std::pair< double, double > xPDF
Definition: PdfInfo.h:14
HepMC::FourVector FourVector(const reco::Candidate::LorentzVector &lvec)
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual const Point & vertex() const =0
vertex position
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int signalProcessID() const
const_reference at(size_type pos) const
std::pair< int, int > id
Definition: PdfInfo.h:12
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Run.h:318
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
const PDF * pdf() const
std::vector< double > & weights()
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
GenParticles2HepMCConverter(const edm::ParameterSet &pset)
double scalePDF
Definition: PdfInfo.h:15
def move(src, dest)
Definition: eostools.py:511
*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
Definition: Run.h:45