CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MergedGenParticleProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
3 #include "MergedGenParticleProducer.hh"
4 
6 
10 
11 #include "HepPDT/ParticleID.hh"
12 
13 MergedGenParticleProducer::MergedGenParticleProducer(const edm::ParameterSet& config) {
14  input_pruned_ = consumes<edm::View<reco::GenParticle>>(config.getParameter<edm::InputTag>("inputPruned"));
15  input_packed_ = consumes<edm::View<pat::PackedGenParticle>>(config.getParameter<edm::InputTag>("inputPacked"));
16 
17  produces<reco::GenParticleCollection>();
18 }
19 
20 void MergedGenParticleProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
21  // Need a ref to the product now for creating the mother/daughter refs
22  auto ref = event.getRefBeforePut<reco::GenParticleCollection>();
23 
24  // Get the input collections
26  event.getByToken(input_pruned_, pruned_handle);
27 
29  event.getByToken(input_packed_, packed_handle);
30 
31  // First determine which packed particles are also still in the pruned collection
32  // so that we can skip them later
33  std::map<pat::PackedGenParticle const*, reco::GenParticle const*> st1_dup_map;
34 
35  // Also map pointers in the original pruned collection to their index in the vector.
36  // This index will be the same in the merged collection.
37  std::map<reco::Candidate const*, std::size_t> pruned_idx_map;
38 
39  for (unsigned int i = 0; i < pruned_handle->size(); ++i) {
40  reco::GenParticle const& src = pruned_handle->at(i);
41  pruned_idx_map[&src] = i;
42  if (src.status() != 1)
43  continue;
44 
45  // Convert the pruned GenParticle into a PackedGenParticle then do an exact
46  // floating point comparison of the pt, eta and phi
47  // This of course relies on the PackedGenParticle constructor being identical
48  // between the CMSSW version this sample was produced with and the one we're
49  // analysing with
51  unsigned found_matches = 0;
52  for (unsigned j = 0; j < packed_handle->size(); ++j) {
53  pat::PackedGenParticle const& pk = packed_handle->at(j);
54  if (pks.pdgId() != pk.pdgId() or pks.p4() != pk.p4())
55  continue;
56  ++found_matches;
57  st1_dup_map[&pk] = &src;
58  }
59  if (found_matches > 1) {
60  edm::LogWarning("MergedGenParticleProducer") << "Found multiple packed matches for: " << i << "\t" << src.pdgId()
61  << "\t" << src.pt() << "\t" << src.y() << "\n";
62  } else if (found_matches == 0 && std::abs(src.y()) < 6.0) {
63  edm::LogWarning("MergedGenParticleProducer")
64  << "unmatched status 1: " << i << "\t" << src.pdgId() << "\t" << src.pt() << "\t" << src.y() << "\n";
65  }
66  }
67 
68  // Fix by Markus
69  // check for photons from pruned (light) hadrons
70  unsigned int nPhotonsFromPrunedHadron = 0;
71  for (unsigned int j = 0; j < packed_handle->size(); ++j) {
72  pat::PackedGenParticle const& pk = packed_handle->at(j);
73  if (isPhotonFromPrunedHadron(pk))
74  ++nPhotonsFromPrunedHadron;
75  }
76 
77  // At this point we know what the size of the merged GenParticle will be so we can create it
78  const unsigned int n =
79  pruned_handle->size() + (packed_handle->size() - st1_dup_map.size()) + nPhotonsFromPrunedHadron;
80  auto cands = std::make_unique<reco::GenParticleCollection>(n);
81 
82  // First copy in all the pruned candidates
83  for (unsigned i = 0; i < pruned_handle->size(); ++i) {
84  reco::GenParticle const& old_cand = pruned_handle->at(i);
85  reco::GenParticle& new_cand = cands->at(i);
86  new_cand = reco::GenParticle(pruned_handle->at(i));
87  // Update the mother and daughter refs to this new merged collection
88  new_cand.resetMothers(ref.id());
89  new_cand.resetDaughters(ref.id());
90  for (unsigned m = 0; m < old_cand.numberOfMothers(); ++m) {
91  new_cand.addMother(reco::GenParticleRef(ref, pruned_idx_map.at(old_cand.mother(m))));
92  }
93  for (unsigned d = 0; d < old_cand.numberOfDaughters(); ++d) {
94  new_cand.addDaughter(reco::GenParticleRef(ref, pruned_idx_map.at(old_cand.daughter(d))));
95  }
96  }
97 
98  // Now copy in the packed candidates that are not already in the pruned
99  for (unsigned i = 0, idx = pruned_handle->size(); i < packed_handle->size(); ++i) {
100  pat::PackedGenParticle const& pk = packed_handle->at(i);
101  if (st1_dup_map.count(&pk))
102  continue;
103  reco::GenParticle& new_cand = cands->at(idx);
104  new_cand = reco::GenParticle(pk.charge(), pk.p4(), pk.vertex(), pk.pdgId(), 1, true);
105 
106  // Insert dummy pi0 mothers for orphaned photons
107  if (isPhotonFromPrunedHadron(pk)) {
108  ++idx;
109  reco::GenParticle& dummy_mother = cands->at(idx);
110  dummy_mother = reco::GenParticle(0, pk.p4(), pk.vertex(), 111, 2, true);
111  for (unsigned m = 0; m < pk.numberOfMothers(); ++m) {
112  new_cand.addMother(reco::GenParticleRef(ref, idx));
113  // Since the packed candidates drop the vertex position we'll take this from the mother
114  if (m == 0) {
115  dummy_mother.setP4(pk.mother(m)->p4());
116  dummy_mother.setVertex(pk.mother(m)->vertex());
117  new_cand.setVertex(pk.mother(m)->vertex());
118  }
119  // Should then add the GenParticle as a daughter of its dummy mother
120  dummy_mother.addDaughter(reco::GenParticleRef(ref, idx - 1));
121  }
122  }
123  // Connect to mother from pruned particles
124  reco::GenParticle& daughter = cands->at(idx);
125  for (unsigned m = 0; m < pk.numberOfMothers(); ++m) {
126  daughter.addMother(reco::GenParticleRef(ref, pruned_idx_map.at(pk.mother(m))));
127  // Since the packed candidates drop the vertex position we'll take this from the mother
128  if (m == 0) {
129  daughter.setVertex(pk.mother(m)->vertex());
130  }
131  // Should then add this GenParticle as a daughter of its mother
132  cands->at(pruned_idx_map.at(pk.mother(m))).addDaughter(reco::GenParticleRef(ref, idx));
133  }
134  ++idx;
135  }
136 
137  event.put(std::move(cands));
138 }
139 
140 bool MergedGenParticleProducer::isPhotonFromPrunedHadron(const pat::PackedGenParticle& pk) const {
141  if (pk.pdgId() == 22 and pk.statusFlags().isDirectHadronDecayProduct()) {
142  // no mother
143  if (pk.numberOfMothers() == 0)
144  return true;
145  // miniaod mother not compatible with the status flag
146  HepPDT::ParticleID motherid(pk.mother(0)->pdgId());
147  if (not(motherid.isHadron() and pk.mother(0)->status() == 2))
148  return true;
149  }
150  return false;
151 }
152 
154 DEFINE_FWK_MODULE(MergedGenParticleProducer);
const Point & vertex() const override
vertex position
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
const Candidate * mother(size_type=0) const override
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
double y() const final
rapidity
double pt() const final
transverse momentum
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::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int charge() const override
electric charge
virtual int status() const =0
status word
size_t numberOfDaughters() const override
number of daughters
int status() const final
status word
size_t numberOfMothers() const override
number of mothers
void setVertex(const Point &vertex) override
set vertex
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
int pdgId() const final
PDG identifier.
tuple d
Definition: ztail.py:151
void resetDaughters(const edm::ProductID &id)
set daughters product ID
const LorentzVector & p4() const override
four-momentum Lorentz vecto r
size_t numberOfMothers() const override
number of mothers
bool isDirectHadronDecayProduct() const
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int pdgId() const override
PDG identifier.
virtual const Point & vertex() const =0
vertex position
void addMother(const typename mothers::value_type &)
add a daughter via a reference
virtual int pdgId() const =0
PDG identifier.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const reco::GenStatusFlags & statusFlags() const
tuple config
parse the configuration file
const reco::Candidate * mother(size_type) const override
return mother at a given position (throws an exception)
Log< level::Warning, false > LogWarning
void resetMothers(const edm::ProductID &id)
set mother product ID
void setP4(const LorentzVector &p4) final
set 4-momentum
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector