CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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") : edm::InputTag("prunedGenParticles")),
52  partonicFinalState(params.getParameter<bool>("partonicFinalState")),
53  excludeResonances(params.getParameter<bool>("excludeResonances")),
54  tausAsJets(params.getParameter<bool>("tausAsJets")),
55  ptMin(0.0){
56  if (params.exists("ignoreParticleIDs"))
57  setIgnoredParticles(params.getParameter<std::vector<unsigned int> >
58  ("ignoreParticleIDs"));
59  setExcludeFromResonancePids(params.getParameter<std::vector<unsigned int> >
60  ("excludeFromResonancePids"));
61  isMiniAOD = ( params.exists("isMiniAOD") ? params.getParameter<bool>("isMiniAOD") : (inTag.label()=="packedGenParticles") );
62 
64  edm::LogError("PartonicFinalStateFromMiniAOD") << "Partonic final state not supported for MiniAOD. Falling back to the stable particle selection.";
65  partonicFinalState = false;
66  }
67 
68  produces <reco::CandidatePtrVector> ();
69 
70  input_genpartcoll_token_ = consumes<reco::CandidateView>(inTag);
71  if(isMiniAOD)
72  input_prunedgenpartcoll_token_ = consumes<reco::CandidateView>(prunedInTag);
73 
74 }
75 
77 
78 void InputGenJetsParticleSelector::setIgnoredParticles(const std::vector<unsigned int> &particleIDs)
79 {
80  ignoreParticleIDs = particleIDs;
81  std::sort(ignoreParticleIDs.begin(), ignoreParticleIDs.end());
82 }
83 
84 void InputGenJetsParticleSelector::setExcludeFromResonancePids(const std::vector<unsigned int> &particleIDs)
85 {
86  excludeFromResonancePids = particleIDs;
87  std::sort( excludeFromResonancePids.begin(), excludeFromResonancePids.end());
88 }
89 
91  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
92  return (pdgId > 0 && pdgId < 6) ||
93  pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21;
94  // tops are not considered "regular" partons
95  // but taus eventually are (since they may hadronize later)
96 }
97 
99 {
100  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
101  return (pdgId > 100 && pdgId < 900) ||
102  (pdgId > 1000 && pdgId < 9000);
103 }
104 
106 {
107  // gauge bosons and tops
108  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
109  return (pdgId > 21 && pdgId <= 42) || pdgId == 6 || pdgId == 7 || pdgId == 8 ; //BUG! was 21. 22=gamma..
110 }
111 
113 {
114  pdgId = pdgId > 0 ? pdgId : -pdgId;
115  std::vector<unsigned int>::const_iterator pos =
116  std::lower_bound(ignoreParticleIDs.begin(),
117  ignoreParticleIDs.end(),
118  (unsigned int)pdgId);
119  return pos != ignoreParticleIDs.end() && *pos == (unsigned int)pdgId;
120 }
121 
123 {
124  pdgId = pdgId > 0 ? pdgId : -pdgId;
125  std::vector<unsigned int>::const_iterator pos =
126  std::lower_bound(excludeFromResonancePids.begin(),
128  (unsigned int)pdgId);
129  return pos != excludeFromResonancePids.end() && *pos == (unsigned int)pdgId;
130 
131 }
132 
134  const reco::Candidate *particle)
135 {
136  InputGenJetsParticleSelector::ParticleVector::const_iterator pos =
137  std::lower_bound(p.begin(), p.end(), particle);
138  if (pos == p.end() || *pos != particle)
139  throw cms::Exception("CorruptedData")
140  << "reco::GenEvent corrupted: Unlisted particles"
141  " in decay tree." << std::endl;
142 
143  return pos - p.begin();
144 }
145 
148  const reco::Candidate *particle)
149 {
150  unsigned int npart=particle->numberOfDaughters();
151  if (!npart) return;
152 
153  for (unsigned int i=0;i<npart;++i){
154  unsigned int idx=partIdx(p,particle->daughter(i));
155  if (invalid[idx])
156  continue;
157  invalid[idx] = true;
158  //cout<<"Invalidated: ["<<setw(4)<<idx<<"] With pt:"<<particle->daughter(i)->pt()<<endl;
159  invalidateTree(invalid, p, particle->daughter(i));
160  }
161 }
162 
163 
167  const reco::Candidate *particle) const
168 {
169  unsigned int npart=particle->numberOfDaughters();
170  if (!npart) {return 0;}
171 
172  for (unsigned int i=0;i<npart;++i){
173  unsigned int idx = partIdx(p,particle->daughter(i));
174  if (invalid[idx])
175  continue;
176  if (isParton((particle->daughter(i)->pdgId()))){
177  return 1;
178  }
179  if (isHadron((particle->daughter(i)->pdgId()))){
180  return -1;
181  }
182  int result = testPartonChildren(invalid,p,particle->daughter(i));
183  if (result) return result;
184  }
185  return 0;
186 }
187 
190  const ParticleVector &p,
191  const reco::Candidate *particle) const
192 {
193  unsigned int idx = partIdx(p, particle);
194  int id = particle->pdgId();
195 
196  if (invalid[idx]) return kIndirect;
197 
198  if (isResonance(id) && (particle->status() == 3 || particle->status() == 22) ){
199  return kDirect;
200  }
201 
202 
203  if (!isIgnored(id) && (isParton(id)))
204  return kNo;
205 
206 
207 
208  unsigned int nMo=particle->numberOfMothers();
209  if (!nMo)
210  return kNo;
211 
212 
213  for(unsigned int i=0;i<nMo;++i){
214  ResonanceState result = fromResonance(invalid,p,particle->mother(i));
215  switch(result) {
216  case kNo:
217  break;
218  case kDirect:
219  if (particle->mother(i)->pdgId()==id || isResonance(id))
220  return kDirect;
221  if(!isExcludedFromResonance(id))
222  break;
223  case kIndirect:
224  return kIndirect;
225  }
226  }
227 return kNo;
228 }
229 
230 
232  const ParticleVector &p,
233  const reco::Candidate *particle) const {
234  return testPartonChildren(invalid, p, particle) > 0;
235 }
236 
237 //######################################################
238 //function NEEDED and called per EVENT by FRAMEWORK:
240 
241 
242  std::auto_ptr<reco::CandidatePtrVector> selected_ (new reco::CandidatePtrVector);
243 
244  ParticleVector particles;
245 
247  if(isMiniAOD)
248  {
249  evt.getByToken(input_prunedgenpartcoll_token_, prunedGenParticles );
250 
251  for (edm::View<reco::Candidate>::const_iterator iter=prunedGenParticles->begin();iter!=prunedGenParticles->end();++iter)
252  {
253  if(iter->status()!=1) // to avoid double-counting, skipping stable particles already contained in the collection of PackedGenParticles
254  particles.push_back(&*iter);
255  }
256  }
257 
259  evt.getByToken(input_genpartcoll_token_, genParticles );
260 
261  std::map<const reco::Candidate*,size_t> particlePtrIdxMap;
262 
263  for (edm::View<reco::Candidate>::const_iterator iter=genParticles->begin();iter!=genParticles->end();++iter){
264  particles.push_back(&*iter);
265  particlePtrIdxMap[&*iter] = (iter - genParticles->begin());
266  }
267 
268  std::sort(particles.begin(), particles.end());
269  unsigned int size = particles.size();
270 
271  ParticleBitmap selected(size, false);
272  ParticleBitmap invalid(size, false);
273 
274  for(unsigned int i = 0; i < size; i++) {
275  const reco::Candidate *particle = particles[i];
276  if (invalid[i])
277  continue;
278  if (particle->status() == 1) // selecting stable particles
279  selected[i] = true;
280  if (partonicFinalState && isParton(particle->pdgId())) {
281 
282  if (particle->numberOfDaughters()==0 &&
283  particle->status() != 1) {
284  // some brokenness in event...
285  invalid[i] = true;
286  }
287  else if (!hasPartonChildren(invalid, particles,
288  particle)) {
289  selected[i] = true;
290  invalidateTree(invalid, particles,particle); //this?!?
291  }
292  }
293 
294  }
295 
296  for(size_t idx = 0; idx < size; ++idx){
297  const reco::Candidate *particle = particles[idx];
298  if (!selected[idx] || invalid[idx]){
299  continue;
300  }
301 
302  if (excludeResonances &&
303  fromResonance(invalid, particles, particle)) {
304  invalid[idx] = true;
305  //cout<<"[INPUTSELECTOR] Invalidates FROM RESONANCE!: ["<<setw(4)<<idx<<"] "<<particle->pdgId()<<" "<<particle->pt()<<endl;
306  continue;
307  }
308 
309  if (isIgnored(particle->pdgId())){
310  continue;
311  }
312 
313 
314  if (particle->pt() >= ptMin){
315  selected_->push_back(genParticles->ptrAt(particlePtrIdxMap[particle]));
316  //cout<<"Finally we have: ["<<setw(4)<<idx<<"] "<<setw(4)<<particle->pdgId()<<" "<<particle->pt()<<endl;
317  }
318  }
319  evt.put(selected_);
320 }
321 
322 
323 
324 //define this as a plug-in
bool hasPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const reco::Candidate *particle) const
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle)
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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
virtual double pt() const =0
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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)
bool exists(std::string const &parameterName) const
checks if a parameter exists
double npart
Definition: HydjetWrapper.h:49
std::vector< unsigned int > excludeFromResonancePids
virtual void produce(edm::Event &evt, const edm::EventSetup &evtSetup)
int testPartonChildren(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_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
tuple result
Definition: query.py:137
ResonanceState fromResonance(ParticleBitmap &invalid, const ParticleVector &p, const reco::Candidate *particle) const
virtual int pdgId() const =0
PDG identifier.
bool isParton(const reco::Candidate &c)
Definition: CandMCTag.cc:48
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
static void invalidateTree(InputGenJetsParticleSelector::ParticleBitmap &invalid, const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle)
std::string const & label() const
Definition: InputTag.h:36
std::vector< unsigned int > ignoreParticleIDs
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
std::vector< const reco::Candidate * > ParticleVector
void setExcludeFromResonancePids(const std::vector< unsigned int > &particleIDs)
tuple size
Write out results.