CMS 3D CMS Logo

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  unsigned int nLeptonsFromPrunedPhoton = 0;
40 
41  for (unsigned int i = 0, idx = 0; i < pruned_handle->size(); ++i) {
42  reco::GenParticle const& src = pruned_handle->at(i);
43  pruned_idx_map[&src] = idx;
44  ++idx;
45 
46  // check for electrons+muons from pruned photons
47  if (isLeptonFromPrunedPhoton(src)) {
48  ++nLeptonsFromPrunedPhoton;
49  ++idx;
50  }
51 
52  if (src.status() != 1)
53  continue;
54 
55  // Convert the pruned GenParticle into a PackedGenParticle then do an exact
56  // floating point comparison of the pt, eta and phi
57  // This of course relies on the PackedGenParticle constructor being identical
58  // between the CMSSW version this sample was produced with and the one we're
59  // analysing with
61  unsigned found_matches = 0;
62  for (unsigned j = 0; j < packed_handle->size(); ++j) {
63  pat::PackedGenParticle const& pk = packed_handle->at(j);
64  if (pks.pdgId() != pk.pdgId() or pks.p4() != pk.p4())
65  continue;
66  ++found_matches;
67  st1_dup_map[&pk] = &src;
68  }
69  if (found_matches > 1) {
70  edm::LogWarning("MergedGenParticleProducer") << "Found multiple packed matches for: " << i << "\t" << src.pdgId()
71  << "\t" << src.pt() << "\t" << src.y() << "\n";
72  } else if (found_matches == 0 && std::abs(src.y()) < 6.0) {
73  edm::LogWarning("MergedGenParticleProducer")
74  << "unmatched status 1: " << i << "\t" << src.pdgId() << "\t" << src.pt() << "\t" << src.y() << "\n";
75  }
76  }
77 
78  // Fix by Markus
79  // check for photons from pruned (light) hadrons
80  unsigned int nPhotonsFromPrunedHadron = 0;
81  for (unsigned int j = 0; j < packed_handle->size(); ++j) {
82  pat::PackedGenParticle const& pk = packed_handle->at(j);
83  if (isPhotonFromPrunedHadron(pk))
84  ++nPhotonsFromPrunedHadron;
85  }
86 
87  // At this point we know what the size of the merged GenParticle will be so we can create it
88  const unsigned int n = pruned_handle->size() + (packed_handle->size() - st1_dup_map.size()) +
89  nPhotonsFromPrunedHadron + nLeptonsFromPrunedPhoton;
90  auto cands = std::make_unique<reco::GenParticleCollection>(n);
91 
92  // First copy in all the pruned candidates
93  unsigned idx = 0;
94  for (unsigned i = 0; i < pruned_handle->size(); ++i) {
95  reco::GenParticle const& old_cand = pruned_handle->at(i);
96  reco::GenParticle& new_cand = cands->at(idx);
97  new_cand = reco::GenParticle(pruned_handle->at(i));
98  // Update the mother and daughter refs to this new merged collection
99  new_cand.resetMothers(ref.id());
100  new_cand.resetDaughters(ref.id());
101  // Insert dummy photon mothers for orphaned electrons+muons
102  if (isLeptonFromPrunedPhoton(old_cand)) {
103  ++idx;
104  reco::GenParticle& dummy_mother = cands->at(idx);
105  dummy_mother = reco::GenParticle(0, old_cand.p4(), old_cand.vertex(), 22, 2, true);
106  for (unsigned m = 0; m < old_cand.numberOfMothers(); ++m) {
107  new_cand.addMother(reco::GenParticleRef(ref, idx));
108  // Since the packed candidates drop the vertex position we'll take this from the mother
109  if (m == 0) {
110  dummy_mother.setP4(old_cand.mother(0)->p4());
111  dummy_mother.setVertex(old_cand.mother(0)->vertex());
112  new_cand.setVertex(old_cand.mother(0)->vertex());
113  }
114  // Should then add the GenParticle as a daughter of its dummy mother
115  dummy_mother.addDaughter(reco::GenParticleRef(ref, idx - 1));
116  for (unsigned m = 0; m < old_cand.numberOfMothers(); ++m) {
117  dummy_mother.addMother(reco::GenParticleRef(ref, pruned_idx_map.at(old_cand.mother(m))));
118  }
119  }
120  } else {
121  for (unsigned m = 0; m < old_cand.numberOfMothers(); ++m) {
122  new_cand.addMother(reco::GenParticleRef(ref, pruned_idx_map.at(old_cand.mother(m))));
123  }
124  }
125  for (unsigned d = 0; d < old_cand.numberOfDaughters(); ++d) {
126  new_cand.addDaughter(reco::GenParticleRef(ref, pruned_idx_map.at(old_cand.daughter(d))));
127  }
128  ++idx;
129  }
130 
131  // Now copy in the packed candidates that are not already in the pruned
132  for (unsigned i = 0; i < packed_handle->size(); ++i) {
133  pat::PackedGenParticle const& pk = packed_handle->at(i);
134  if (st1_dup_map.count(&pk))
135  continue;
136  reco::GenParticle& new_cand = cands->at(idx);
137  new_cand = reco::GenParticle(pk.charge(), pk.p4(), pk.vertex(), pk.pdgId(), 1, true);
138 
139  // Insert dummy pi0 mothers for orphaned photons
140  if (isPhotonFromPrunedHadron(pk)) {
141  ++idx;
142  reco::GenParticle& dummy_mother = cands->at(idx);
143  dummy_mother = reco::GenParticle(0, pk.p4(), pk.vertex(), 111, 2, true);
144  for (unsigned m = 0; m < pk.numberOfMothers(); ++m) {
145  new_cand.addMother(reco::GenParticleRef(ref, idx));
146  // Since the packed candidates drop the vertex position we'll take this from the mother
147  if (m == 0) {
148  dummy_mother.setP4(pk.mother(m)->p4());
149  dummy_mother.setVertex(pk.mother(m)->vertex());
150  new_cand.setVertex(pk.mother(m)->vertex());
151  }
152  // Should then add the GenParticle as a daughter of its dummy mother
153  dummy_mother.addDaughter(reco::GenParticleRef(ref, idx - 1));
154  }
155  }
156  // Connect to mother from pruned particles
157  reco::GenParticle& daughter = cands->at(idx);
158  for (unsigned m = 0; m < pk.numberOfMothers(); ++m) {
159  daughter.addMother(reco::GenParticleRef(ref, pruned_idx_map.at(pk.mother(m))));
160  // Since the packed candidates drop the vertex position we'll take this from the mother
161  if (m == 0) {
162  daughter.setVertex(pk.mother(m)->vertex());
163  }
164  // Should then add this GenParticle as a daughter of its mother
165  cands->at(pruned_idx_map.at(pk.mother(m))).addDaughter(reco::GenParticleRef(ref, idx));
166  }
167  ++idx;
168  }
169 
170  event.put(std::move(cands));
171 }
172 
173 bool MergedGenParticleProducer::isPhotonFromPrunedHadron(const pat::PackedGenParticle& pk) const {
174  if (pk.pdgId() == 22 and pk.statusFlags().isDirectHadronDecayProduct()) {
175  // no mother
176  if (pk.numberOfMothers() == 0)
177  return true;
178  // miniaod mother not compatible with the status flag
179  HepPDT::ParticleID motherid(pk.mother(0)->pdgId());
180  if (not(motherid.isHadron() and pk.mother(0)->status() == 2))
181  return true;
182  }
183  return false;
184 }
185 
186 bool MergedGenParticleProducer::isLeptonFromPrunedPhoton(const reco::GenParticle& pk) const {
187  if ((abs(pk.pdgId()) == 11 or abs(pk.pdgId()) == 13) and
189  // this is probably not a prompt lepton but from pair production via a pruned photon
190  if (pk.numberOfMothers() > 0 and pk.mother(0)->pdgId() != 22) {
191  return true;
192  }
193  }
194  return false;
195 }
196 
198 DEFINE_FWK_MODULE(MergedGenParticleProducer);
bool fromHardProcess() const
const Point & vertex() const override
vertex position
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool isDirectHadronDecayProduct() const
const Candidate * mother(size_type=0) const override
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
int charge() const override
electric charge
virtual int status() const =0
status word
size_t numberOfDaughters() const override
number of daughters
const Point & vertex() const override
vertex position (overwritten by PF...)
Definition: config.py:1
size_t numberOfMothers() const override
number of mothers
const LorentzVector & p4() const final
four-momentum Lorentz vector
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.
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
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
int pdgId() const override
PDG identifier.
virtual const Point & vertex() const =0
vertex position
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
d
Definition: ztail.py:151
void addMother(const typename mothers::value_type &)
add a daughter via a reference
virtual int pdgId() const =0
PDG identifier.
const GenStatusFlags & statusFlags() const
Definition: GenParticle.h:38
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
bool isDirectTauDecayProduct() const
def move(src, dest)
Definition: eostools.py:511
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
Definition: event.py:1
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
const reco::GenStatusFlags & statusFlags() const