CMS 3D CMS Logo

GenParticleProducer.cc
Go to the documentation of this file.
1 /* \class GenParticleProducer
2  *
3  * \author Luca Lista, INFN
4  *
5  * Convert HepMC GenEvent format into a collection of type
6  * CandidateCollection containing objects of type GenParticle
7  *
8  *
9  */
20 
21 #include <vector>
22 #include <string>
23 #include <unordered_map>
24 
25 namespace edm {
26  class ParameterSet;
27 }
28 namespace HepMC {
29  class GenParticle;
30  class GenEvent;
31 } // namespace HepMC
32 
33 static constexpr int PDGCacheMax = 32768;
34 
35 namespace {
36  struct IDto3Charge {
37  IDto3Charge(HepPDT::ParticleDataTable const&, bool abortOnUnknownPDGCode);
38 
39  int chargeTimesThree(int) const;
40 
41  private:
42  std::vector<int> chargeP_, chargeM_;
43  std::unordered_map<int, int> chargeMap_;
44  bool abortOnUnknownPDGCode_;
45  };
46 
47  IDto3Charge::IDto3Charge(HepPDT::ParticleDataTable const& iTable, bool iAbortOnUnknownPDGCode)
48  : chargeP_(PDGCacheMax, 0), chargeM_(PDGCacheMax, 0), abortOnUnknownPDGCode_(iAbortOnUnknownPDGCode) {
49  for (auto const& p : iTable) {
50  const HepPDT::ParticleID& id = p.first;
51  int pdgId = id.pid(), apdgId = std::abs(pdgId);
52  int q3 = id.threeCharge();
53  if (apdgId < PDGCacheMax && pdgId > 0) {
54  chargeP_[apdgId] = q3;
55  chargeM_[apdgId] = -q3;
56  } else if (apdgId < PDGCacheMax) {
57  chargeP_[apdgId] = -q3;
58  chargeM_[apdgId] = q3;
59  } else {
60  chargeMap_.emplace(pdgId, q3);
61  chargeMap_.emplace(-pdgId, -q3);
62  }
63  }
64  }
65 
66  int IDto3Charge::chargeTimesThree(int id) const {
67  if (std::abs(id) < PDGCacheMax)
68  return id > 0 ? chargeP_[id] : chargeM_[-id];
69  auto f = chargeMap_.find(id);
70  if (f == chargeMap_.end()) {
71  if (abortOnUnknownPDGCode_)
72  throw edm::Exception(edm::errors::LogicError) << "invalid PDG id: " << id;
73  else
74  return HepPDT::ParticleID(id).threeCharge();
75  }
76  return f->second;
77  }
78 
79 } // namespace
80 
81 class GenParticleProducer : public edm::global::EDProducer<edm::RunCache<IDto3Charge>> {
82 public:
86  ~GenParticleProducer() override;
87 
89  void produce(edm::StreamID, edm::Event& e, const edm::EventSetup&) const override;
90  std::shared_ptr<IDto3Charge> globalBeginRun(const edm::Run&, const edm::EventSetup&) const override;
91  void globalEndRun(edm::Run const&, edm::EventSetup const&) const override{};
92 
93  bool convertParticle(reco::GenParticle& cand, const HepMC::GenParticle* part, const IDto3Charge& id2Charge) const;
95  const HepMC::GenParticle* part,
96  reco::GenParticleRefProd const& ref,
97  size_t index,
98  std::unordered_map<int, size_t>& barcodes) const;
99  bool fillIndices(const HepMC::GenEvent* mc,
100  std::vector<const HepMC::GenParticle*>& particles,
101  std::vector<int>& barCodeVector,
102  int offset,
103  std::unordered_map<int, size_t>& barcodes) const;
104 
105 private:
108  std::vector<edm::EDGetTokenT<edm::HepMCProduct>> vectorSrcTokens_;
111 
116 
119  bool useCF_;
120 
123 };
124 
126 
127 //#include "SimDataFormats/HiGenData/interface/SubEventMap.h"
129 
138 
139 #include <fstream>
140 #include <algorithm>
141 using namespace edm;
142 using namespace reco;
143 using namespace std;
144 using namespace HepMC;
145 
146 static const double mmToCm = 0.1;
147 static const double mmToNs = 1.0 / 299792458e-6;
148 
150  : abortOnUnknownPDGCode_(false),
151  saveBarCodes_(cfg.getUntrackedParameter<bool>("saveBarCodes", false)),
152  doSubEvent_(cfg.getUntrackedParameter<bool>("doSubEvent", false)),
153  useCF_(cfg.getUntrackedParameter<bool>("useCrossingFrame", false)) {
154  particleTableToken_ = esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord, edm::Transition::BeginRun>();
155  produces<GenParticleCollection>();
156  produces<math::XYZPointF>("xyz0");
157  produces<float>("t0");
158  if (saveBarCodes_) {
159  std::string alias(cfg.getParameter<std::string>("@module_label"));
160  produces<vector<int>>().setBranchAlias(alias + "BarCodes");
161  }
162 
163  if (useCF_)
164  mixToken_ =
165  mayConsume<CrossingFrame<HepMCProduct>>(InputTag(cfg.getParameter<std::string>("mix"), "generatorSmeared"));
166  else
167  srcToken_ = mayConsume<HepMCProduct>(cfg.getParameter<InputTag>("src"));
168 }
169 
171 
172 std::shared_ptr<IDto3Charge> GenParticleProducer::globalBeginRun(const Run&, const EventSetup& es) const {
174  return std::make_shared<IDto3Charge>(*pdt, abortOnUnknownPDGCode_);
175 }
176 
178  std::unordered_map<int, size_t> barcodes;
179 
180  size_t totalSize = 0;
181  const GenEvent* mc = nullptr;
182  MixCollection<HepMCProduct>* cfhepmcprod = nullptr;
183  size_t npiles = 1;
184 
185  if (useCF_) {
187  evt.getByToken(mixToken_, cf);
188  cfhepmcprod = new MixCollection<HepMCProduct>(cf.product());
189  npiles = cfhepmcprod->size();
190  LogDebug("GenParticleProducer") << "npiles : " << npiles << endl;
191  for (unsigned int icf = 0; icf < npiles; ++icf) {
192  LogDebug("GenParticleProducer") << "subSize : " << cfhepmcprod->getObject(icf).GetEvent()->particles_size()
193  << endl;
194  totalSize += cfhepmcprod->getObject(icf).GetEvent()->particles_size();
195  }
196  LogDebug("GenParticleProducer") << "totalSize : " << totalSize << endl;
197  } else {
199  evt.getByToken(srcToken_, mcp);
200  mc = mcp->GetEvent();
201  if (mc == nullptr)
202  throw edm::Exception(edm::errors::InvalidReference) << "HepMC has null pointer to GenEvent" << endl;
203  totalSize = mc->particles_size();
204  }
205 
206  // initialise containers
207  const size_t size = totalSize;
208  vector<const HepMC::GenParticle*> particles(size);
209  auto candsPtr = std::make_unique<GenParticleCollection>(size);
210  auto barCodeVector = std::make_unique<vector<int>>(size);
211  std::unique_ptr<math::XYZPointF> xyz0Ptr(new math::XYZPointF(0., 0., 0.));
212  std::unique_ptr<float> t0Ptr(new float(0.f));
214  GenParticleCollection& cands = *candsPtr;
215  size_t offset = 0;
216  size_t suboffset = 0;
217 
218  IDto3Charge const& id2Charge = *runCache(evt.getRun().index());
220  if (doSubEvent_ || useCF_) {
221  for (size_t ipile = 0; ipile < npiles; ++ipile) {
222  LogDebug("GenParticleProducer") << "mixed object ipile : " << ipile << endl;
223  barcodes.clear();
224  if (useCF_)
225  mc = cfhepmcprod->getObject(ipile).GetEvent();
226 
227  //Look whether heavy ion/signal event
228  bool isHI = false;
229  const HepMC::HeavyIon* hi = mc->heavy_ion();
230  if (hi && hi->Ncoll_hard() > 1)
231  isHI = true;
232  size_t num_particles = mc->particles_size();
233  LogDebug("GenParticleProducer") << "num_particles : " << num_particles << endl;
234  if (ipile == 0) {
235  auto origin = (*mc->vertices_begin())->position();
236  xyz0Ptr->SetXYZ(origin.x() * mmToCm, origin.y() * mmToCm, origin.z() * mmToCm);
237  *t0Ptr = origin.t() * mmToNs;
238  }
239  fillIndices(mc, particles, *barCodeVector, offset, barcodes);
240  // fill output collection and save association
241  for (size_t ipar = offset; ipar < offset + num_particles; ++ipar) {
242  const HepMC::GenParticle* part = particles[ipar];
243  reco::GenParticle& cand = cands[ipar];
244  // convert HepMC::GenParticle to new reco::GenParticle
245  convertParticle(cand, part, id2Charge);
246  cand.resetDaughters(ref.id());
247  }
248 
249  for (size_t d = offset; d < offset + num_particles; ++d) {
251  const GenVertex* productionVertex = part->production_vertex();
252  int sub_id = 0;
253  if (productionVertex != nullptr) {
254  sub_id = productionVertex->id();
255  if (!isHI)
256  sub_id = 0;
257  // search barcode map and attach daughters
258  fillDaughters(cands, part, ref, d, barcodes);
259  } else {
260  const GenVertex* endVertex = part->end_vertex();
261  if (endVertex != nullptr)
262  sub_id = endVertex->id();
263  else
264  throw cms::Exception("SubEventID")
265  << "SubEvent not determined. Particle has no production and no end vertex!" << endl;
266  }
267  if (sub_id < 0)
268  sub_id = 0;
269  int new_id = sub_id + suboffset;
270  GenParticleRef dref(ref, d);
271  cands[d].setCollisionId(new_id); // For new GenParticle
272  LogDebug("VertexId") << "SubEvent offset 3 : " << suboffset;
273  }
274  int nsub = -2;
275  if (isHI) {
276  nsub = hi->Ncoll_hard() + 1;
277  suboffset += nsub;
278  } else {
279  suboffset += 1;
280  }
281  offset += num_particles;
282  }
283  } else {
284  auto origin = (*mc->vertices_begin())->position();
285  xyz0Ptr->SetXYZ(origin.x() * mmToCm, origin.y() * mmToCm, origin.z() * mmToCm);
286  *t0Ptr = origin.t() * mmToNs;
287  fillIndices(mc, particles, *barCodeVector, 0, barcodes);
288 
289  // fill output collection and save association
290  for (size_t i = 0; i < particles.size(); ++i) {
293  // convert HepMC::GenParticle to new reco::GenParticle
294  convertParticle(cand, part, id2Charge);
295  cand.resetDaughters(ref.id());
296  }
297 
298  // fill references to daughters
299  for (size_t d = 0; d < cands.size(); ++d) {
301  const GenVertex* productionVertex = part->production_vertex();
302  // search barcode map and attach daughters
303  if (productionVertex != nullptr)
304  fillDaughters(cands, part, ref, d, barcodes);
305  cands[d].setCollisionId(0);
306  }
307  }
308 
309  evt.put(std::move(candsPtr));
310  if (saveBarCodes_)
311  evt.put(std::move(barCodeVector));
312  if (cfhepmcprod)
313  delete cfhepmcprod;
314  evt.put(std::move(xyz0Ptr), "xyz0");
315  evt.put(std::move(t0Ptr), "t0");
316 }
317 
319  const HepMC::GenParticle* part,
320  IDto3Charge const& id2Charge) const {
321  Candidate::LorentzVector p4(part->momentum());
322  int pdgId = part->pdg_id();
323  cand.setThreeCharge(id2Charge.chargeTimesThree(pdgId));
324  cand.setPdgId(pdgId);
325  cand.setStatus(part->status());
326  cand.setP4(p4);
327  cand.setCollisionId(0);
328  const GenVertex* v = part->production_vertex();
329  if (v != nullptr) {
330  ThreeVector vtx = v->point3d();
331  Candidate::Point vertex(vtx.x() * mmToCm, vtx.y() * mmToCm, vtx.z() * mmToCm);
332  cand.setVertex(vertex);
333  } else {
334  cand.setVertex(Candidate::Point(0, 0, 0));
335  }
336  mcTruthHelper_.fillGenStatusFlags(*part, cand.statusFlags());
337  return true;
338 }
339 
341  const HepMC::GenParticle* part,
342  reco::GenParticleRefProd const& ref,
343  size_t index,
344  std::unordered_map<int, size_t>& barcodes) const {
345  const GenVertex* productionVertex = part->production_vertex();
346  size_t numberOfMothers = productionVertex->particles_in_size();
347  if (numberOfMothers > 0) {
348  GenVertex::particles_in_const_iterator motherIt = productionVertex->particles_in_const_begin();
349  for (; motherIt != productionVertex->particles_in_const_end(); motherIt++) {
350  const HepMC::GenParticle* mother = *motherIt;
351  size_t m = barcodes.find(mother->barcode())->second;
352  cands[m].addDaughter(GenParticleRef(ref, index));
353  cands[index].addMother(GenParticleRef(ref, m));
354  }
355  }
356 
357  return true;
358 }
359 
361  vector<const HepMC::GenParticle*>& particles,
362  vector<int>& barCodeVector,
363  int offset,
364  std::unordered_map<int, size_t>& barcodes) const {
365  size_t idx = offset;
366  HepMC::GenEvent::particle_const_iterator begin = mc->particles_begin(), end = mc->particles_end();
367  for (HepMC::GenEvent::particle_const_iterator p = begin; p != end; ++p) {
368  const HepMC::GenParticle* particle = *p;
369  size_t barCode_this_event = particle->barcode();
370  size_t barCode = barCode_this_event + offset;
371  if (barcodes.find(barCode) != barcodes.end())
372  throw cms::Exception("WrongReference") << "barcodes are duplicated! " << endl;
373  particles[idx] = particle;
374  barCodeVector[idx] = barCode;
375  barcodes.insert(make_pair(barCode_this_event, idx++));
376  }
377  return true;
378 }
379 
381 
size
Write out results.
GenParticleProducer(const edm::ParameterSet &)
constructor
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool abortOnUnknownPDGCode_
unknown code treatment flag
MCTruthHelper< HepMC::GenParticle > mcTruthHelper_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
static const double mmToNs
MCTruthHelper< reco::GenParticle > mcTruthHelperGenParts_
HepPDT::ParticleDataTable ParticleDataTable
const T & getObject(unsigned int ip) const
Definition: MixCollection.h:26
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particleTableToken_
T const * product() const
Definition: Handle.h:70
std::vector< edm::EDGetTokenT< edm::HepMCProduct > > vectorSrcTokens_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:536
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > mixToken_
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
edm::EDGetTokenT< edm::HepMCProduct > srcToken_
source collection name
bool fillIndices(const HepMC::GenEvent *mc, std::vector< const HepMC::GenParticle *> &particles, std::vector< int > &barCodeVector, int offset, std::unordered_map< int, size_t > &barcodes) const
const int nsub
U second(std::pair< T, U > const &p)
int size() const
Definition: MixCollection.h:19
Definition: EPCuts.h:4
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Run const & getRun() const
Definition: Event.cc:112
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
d
Definition: ztail.py:151
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
RefProd< PROD > getRefBeforePut()
Definition: Event.h:158
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:124
void globalEndRun(edm::Run const &, edm::EventSetup const &) const override
static constexpr int PDGCacheMax
static const double mmToCm
part
Definition: HCALResponse.h:20
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
~GenParticleProducer() override
destructor
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
void produce(edm::StreamID, edm::Event &e, const edm::EventSetup &) const override
process one event
bool doSubEvent_
input & output modes
RunIndex index() const
Definition: Run.cc:28
fixed size matrix
HLT enums.
std::shared_ptr< IDto3Charge > globalBeginRun(const edm::Run &, const edm::EventSetup &) const override
bool fillDaughters(reco::GenParticleCollection &cand, const HepMC::GenParticle *part, reco::GenParticleRefProd const &ref, size_t index, std::unordered_map< int, size_t > &barcodes) const
static int position[264][3]
Definition: ReadPGInfo.cc:289
void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags) const
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
bool saveBarCodes_
save bar-codes
def move(src, dest)
Definition: eostools.py:511
bool convertParticle(reco::GenParticle &cand, const HepMC::GenParticle *part, const IDto3Charge &id2Charge) const
Definition: Run.h:45
#define LogDebug(id)