CMS 3D CMS Logo

InputGenJetsParticleSelector.cc
Go to the documentation of this file.
1 /* \class GenJetInputParticleSelector
2 *
3 * Selects particles that are used as input for the GenJet collection.
4 * Logic: select all stable particles, except for particles specified in
5 * the config file that come from
6 * W,Z and H decays, and except for a special list, which can be used for
7 * unvisible BSM-particles.
8 * It is also possible to only selected the partonic final state,
9 * which means all particles before the hadronization step.
10 *
11 * The algorithm is based on code of Christophe Saout.
12 *
13 * Usage: [example for no resonance from nu an mu, and deselect invisible BSM
14 * particles ]
15 *
16 * module genJetParticles = InputGenJetsParticleSelector {
17 * InputTag src = "genParticles"
18 * bool partonicFinalState = false
19 * bool excludeResonances = true
20 * vuint32 excludeFromResonancePids = {13,12,14,16}
21 * bool tausAsJets = false
22 * vuint32 ignoreParticleIDs = { 1000022, 2000012, 2000014,
23 * 2000016, 1000039, 5000039,
24 * 4000012, 9900012, 9900014,
25 * 9900016, 39}
26 * }
27 *
28 *
29 * \author: Christophe Saout, Andreas Oehler
30 *
31 * Modifications:
32 *
33 * 04.08.2014: Dinko Ferencek
34 * Added support for Pythia8 (status=22 for intermediate resonances)
35 * 23.09.2014: Dinko Ferencek
36 * Generalized code to work with miniAOD (except for the partonicFinalState which requires AOD)
37 *
38 */
39 
43 //#include <iostream>
44 #include <memory>
46 
47 using namespace std;
48 
50  : inTag(params.getParameter<edm::InputTag>("src")),
51  prunedInTag(params.exists("prunedGenParticles") ? params.getParameter<edm::InputTag>("prunedGenParticles")
52  : edm::InputTag("prunedGenParticles")),
53  partonicFinalState(params.getParameter<bool>("partonicFinalState")),
54  excludeResonances(params.getParameter<bool>("excludeResonances")),
55  tausAsJets(params.getParameter<bool>("tausAsJets")),
56  ptMin(0.0) {
57  if (params.exists("ignoreParticleIDs"))
58  setIgnoredParticles(params.getParameter<std::vector<unsigned int> >("ignoreParticleIDs"));
59  setExcludeFromResonancePids(params.getParameter<std::vector<unsigned int> >("excludeFromResonancePids"));
60  isMiniAOD =
61  (params.exists("isMiniAOD") ? params.getParameter<bool>("isMiniAOD") : (inTag.label() == "packedGenParticles"));
62 
64  edm::LogError("PartonicFinalStateFromMiniAOD")
65  << "Partonic final state not supported for MiniAOD. Falling back to the stable particle selection.";
66  partonicFinalState = false;
67  }
68 
69  produces<reco::CandidatePtrVector>();
70 
71  input_genpartcoll_token_ = consumes<reco::CandidateView>(inTag);
72  if (isMiniAOD)
73  input_prunedgenpartcoll_token_ = consumes<reco::CandidateView>(prunedInTag);
74 }
75 
77 
78 void InputGenJetsParticleSelector::setIgnoredParticles(const std::vector<unsigned int> &particleIDs) {
79  ignoreParticleIDs = particleIDs;
81 }
82 
83 void InputGenJetsParticleSelector::setExcludeFromResonancePids(const std::vector<unsigned int> &particleIDs) {
84  excludeFromResonancePids = particleIDs;
86 }
87 
89  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
90  return (pdgId > 0 && pdgId < 6) || pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21;
91  // tops are not considered "regular" partons
92  // but taus eventually are (since they may hadronize later)
93 }
94 
96  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
97  return (pdgId > 100 && pdgId < 900) || (pdgId > 1000 && pdgId < 9000);
98 }
99 
101  // gauge bosons and tops
102  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
103  return (pdgId > 21 && pdgId <= 42) || pdgId == 6 || pdgId == 7 || pdgId == 8; //BUG! was 21. 22=gamma..
104 }
105 
107  pdgId = pdgId > 0 ? pdgId : -pdgId;
108  std::vector<unsigned int>::const_iterator pos =
109  std::lower_bound(ignoreParticleIDs.begin(), ignoreParticleIDs.end(), (unsigned int)pdgId);
110  return pos != ignoreParticleIDs.end() && *pos == (unsigned int)pdgId;
111 }
112 
114  pdgId = pdgId > 0 ? pdgId : -pdgId;
115  std::vector<unsigned int>::const_iterator pos =
117  return pos != excludeFromResonancePids.end() && *pos == (unsigned int)pdgId;
118 }
119 
120 static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle) {
121  InputGenJetsParticleSelector::ParticleVector::const_iterator pos = std::lower_bound(p.begin(), p.end(), particle);
122  if (pos == p.end() || *pos != particle)
123  throw cms::Exception("CorruptedData") << "reco::GenEvent corrupted: Unlisted particles"
124  " in decay tree."
125  << std::endl;
126 
127  return pos - p.begin();
128 }
129 
132  const reco::Candidate *particle) {
133  unsigned int npart = particle->numberOfDaughters();
134  if (!npart)
135  return;
136 
137  for (unsigned int i = 0; i < npart; ++i) {
138  unsigned int idx = partIdx(p, particle->daughter(i));
139  if (invalid[idx])
140  continue;
141  invalid[idx] = true;
142  //cout<<"Invalidated: ["<<setw(4)<<idx<<"] With pt:"<<particle->daughter(i)->pt()<<endl;
143  invalidateTree(invalid, p, particle->daughter(i));
144  }
145 }
146 
149  const reco::Candidate *particle) const {
150  unsigned int npart = particle->numberOfDaughters();
151  if (!npart) {
152  return 0;
153  }
154 
155  for (unsigned int i = 0; i < npart; ++i) {
156  unsigned int idx = partIdx(p, particle->daughter(i));
157  if (invalid[idx])
158  continue;
159  if (isParton((particle->daughter(i)->pdgId()))) {
160  return 1;
161  }
162  if (isHadron((particle->daughter(i)->pdgId()))) {
163  return -1;
164  }
165  int result = testPartonChildren(invalid, p, particle->daughter(i));
166  if (result)
167  return result;
168  }
169  return 0;
170 }
171 
173  ParticleBitmap &invalid, const ParticleVector &p, const reco::Candidate *particle) const {
174  unsigned int idx = partIdx(p, particle);
175  int id = particle->pdgId();
176 
177  if (invalid[idx])
178  return kIndirect;
179 
180  if (isResonance(id) && (particle->status() == 3 || particle->status() == 22)) {
181  return kDirect;
182  }
183 
184  if (!isIgnored(id) && (isParton(id)))
185  return kNo;
186 
187  unsigned int nMo = particle->numberOfMothers();
188  if (!nMo)
189  return kNo;
190 
191  for (unsigned int i = 0; i < nMo; ++i) {
193  switch (result) {
194  case kNo:
195  break;
196  case kDirect:
197  if (particle->mother(i)->pdgId() == id || isResonance(id))
198  return kDirect;
199  if (!isExcludedFromResonance(id))
200  break;
201  [[fallthrough]];
202  case kIndirect:
203  return kIndirect;
204  }
205  }
206  return kNo;
207 }
208 
210  const ParticleVector &p,
211  const reco::Candidate *particle) const {
212  return testPartonChildren(invalid, p, particle) > 0;
213 }
214 
215 //######################################################
216 //function NEEDED and called per EVENT by FRAMEWORK:
218  auto selected_ = std::make_unique<reco::CandidatePtrVector>();
219 
221 
223  if (isMiniAOD) {
225 
227  iter != prunedGenParticles->end();
228  ++iter) {
229  if (iter->status() !=
230  1) // to avoid double-counting, skipping stable particles already contained in the collection of PackedGenParticles
231  particles.push_back(&*iter);
232  }
233  }
234 
237 
238  std::map<const reco::Candidate *, size_t> particlePtrIdxMap;
239 
240  for (edm::View<reco::Candidate>::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
241  particles.push_back(&*iter);
242  particlePtrIdxMap[&*iter] = (iter - genParticles->begin());
243  }
244 
245  std::sort(particles.begin(), particles.end());
246  unsigned int size = particles.size();
247 
248  ParticleBitmap selected(size, false);
249  ParticleBitmap invalid(size, false);
250 
251  for (unsigned int i = 0; i < size; i++) {
252  const reco::Candidate *particle = particles[i];
253  if (invalid[i])
254  continue;
255  if (particle->status() == 1) // selecting stable particles
256  selected[i] = true;
257  if (partonicFinalState && isParton(particle->pdgId())) {
258  if (particle->numberOfDaughters() == 0 && particle->status() != 1) {
259  // some brokenness in event...
260  invalid[i] = true;
261  } else if (!hasPartonChildren(invalid, particles, particle)) {
262  selected[i] = true;
263  invalidateTree(invalid, particles, particle); //this?!?
264  }
265  }
266  }
267 
268  for (size_t idx = 0; idx < size; ++idx) {
269  const reco::Candidate *particle = particles[idx];
270  if (!selected[idx] || invalid[idx]) {
271  continue;
272  }
273 
274  if (excludeResonances && fromResonance(invalid, particles, particle)) {
275  invalid[idx] = true;
276  //cout<<"[INPUTSELECTOR] Invalidates FROM RESONANCE!: ["<<setw(4)<<idx<<"] "<<particle->pdgId()<<" "<<particle->pt()<<endl;
277  continue;
278  }
279 
280  if (isIgnored(particle->pdgId())) {
281  continue;
282  }
283 
284  if (particle->pt() >= ptMin) {
285  selected_->push_back(genParticles->ptrAt(particlePtrIdxMap[particle]));
286  //cout<<"Finally we have: ["<<setw(4)<<idx<<"] "<<setw(4)<<particle->pdgId()<<" "<<particle->pt()<<endl;
287  }
288  }
289  evt.put(std::move(selected_));
290 }
291 
292 //define this as a plug-in
size
Write out results.
static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
edm::EDGetTokenT< reco::CandidateView > input_genpartcoll_token_
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
virtual double pt() const =0
transverse momentum
virtual int status() const =0
status word
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
constexpr float ptMin
double npart
Definition: HydjetWrapper.h:48
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
std::vector< unsigned int > excludeFromResonancePids
std::string const & label() const
Definition: InputTag.h:36
Log< level::Error, false > LogError
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &evtSetup) const override
bool hasPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const reco::Candidate *particle) const
void setIgnoredParticles(const std::vector< unsigned int > &particleIDs)
virtual size_type numberOfDaughters() const =0
number of daughters
edm::EDGetTokenT< reco::CandidateView > input_prunedgenpartcoll_token_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual int pdgId() const =0
PDG identifier.
int testPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const reco::Candidate *particle) const
static void invalidateTree(InputGenJetsParticleSelector::ParticleBitmap &invalid, const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle)
std::vector< unsigned int > ignoreParticleIDs
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
std::vector< const reco::Candidate * > ParticleVector
void setExcludeFromResonancePids(const std::vector< unsigned int > &particleIDs)
ResonanceState fromResonance(ParticleBitmap &invalid, const ParticleVector &p, const reco::Candidate *particle) const
def move(src, dest)
Definition: eostools.py:511