CMS 3D CMS Logo

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