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 *
36 */
37 
41 //#include <iostream>
42 #include <memory>
44 
45 using namespace std;
46 
48  inTag(params.getParameter<edm::InputTag>("src")),
49  partonicFinalState(params.getParameter<bool>("partonicFinalState")),
50  excludeResonances(params.getParameter<bool>("excludeResonances")),
51  tausAsJets(params.getParameter<bool>("tausAsJets")),
52  ptMin(0.0){
53  if (params.exists("ignoreParticleIDs"))
54  setIgnoredParticles(params.getParameter<std::vector<unsigned int> >
55  ("ignoreParticleIDs"));
56  setExcludeFromResonancePids(params.getParameter<std::vector<unsigned int> >
57  ("excludeFromResonancePids"));
58 
59  produces <reco::GenParticleRefVector> ();
60 
61  input_genpartcoll_token_ = consumes<reco::GenParticleCollection>(inTag);
62 
63 }
64 
66 
67 void InputGenJetsParticleSelector::setIgnoredParticles(const std::vector<unsigned int> &particleIDs)
68 {
69  ignoreParticleIDs = particleIDs;
71 }
72 
73 void InputGenJetsParticleSelector::setExcludeFromResonancePids(const std::vector<unsigned int> &particleIDs)
74 {
75  excludeFromResonancePids = particleIDs;
77 }
78 
80  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
81  return (pdgId > 0 && pdgId < 6) ||
82  pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21;
83  // tops are not considered "regular" partons
84  // but taus eventually are (since they may hadronize later)
85 }
86 
88 {
89  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
90  return (pdgId > 100 && pdgId < 900) ||
91  (pdgId > 1000 && pdgId < 9000);
92 }
93 
95 {
96  // gauge bosons and tops
97  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
98  return (pdgId > 21 && pdgId <= 42) || pdgId == 6 || pdgId == 7 || pdgId == 8 ; //BUG! was 21. 22=gamma..
99 }
100 
102 {
103  pdgId = pdgId > 0 ? pdgId : -pdgId;
104  std::vector<unsigned int>::const_iterator pos =
105  std::lower_bound(ignoreParticleIDs.begin(),
106  ignoreParticleIDs.end(),
107  (unsigned int)pdgId);
108  return pos != ignoreParticleIDs.end() && *pos == (unsigned int)pdgId;
109 }
110 
112 {
113  pdgId = pdgId > 0 ? pdgId : -pdgId;
114  std::vector<unsigned int>::const_iterator pos =
115  std::lower_bound(excludeFromResonancePids.begin(),
117  (unsigned int)pdgId);
118  return pos != excludeFromResonancePids.end() && *pos == (unsigned int)pdgId;
119 
120 }
121 
123  const reco::GenParticle *particle)
124 {
125  InputGenJetsParticleSelector::ParticleVector::const_iterator pos =
126  std::lower_bound(p.begin(), p.end(), particle);
127  if (pos == p.end() || *pos != particle)
128  throw cms::Exception("CorruptedData")
129  << "reco::GenEvent corrupted: Unlisted particles"
130  " in decay tree." << std::endl;
131 
132  return pos - p.begin();
133 }
134 
137  const reco::GenParticle *particle)
138 {
139  unsigned int npart=particle->numberOfDaughters();
140  if (!npart) return;
141 
142  for (unsigned int i=0;i<npart;++i){
143  unsigned int idx=partIdx(p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
144  if (invalid[idx])
145  continue;
146  invalid[idx] = true;
147  //cout<<"Invalidated: ["<<setw(4)<<idx<<"] With pt:"<<particle->daughter(i)->pt()<<endl;
148  invalidateTree(invalid, p, dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
149  }
150 }
151 
152 
156  const reco::GenParticle *particle) const
157 {
158  unsigned int npart=particle->numberOfDaughters();
159  if (!npart) {return 0;}
160 
161  for (unsigned int i=0;i<npart;++i){
162  unsigned int idx = partIdx(p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
163  if (invalid[idx])
164  continue;
165  if (isParton((particle->daughter(i)->pdgId()))){
166  return 1;
167  }
168  if (isHadron((particle->daughter(i)->pdgId()))){
169  return -1;
170  }
171  int result = testPartonChildren(invalid,p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
172  if (result) return result;
173  }
174  return 0;
175 }
176 
179  const ParticleVector &p,
180  const reco::GenParticle *particle) const
181 {
182  unsigned int idx = partIdx(p, particle);
183  int id = particle->pdgId();
184 
185  if (invalid[idx]) return kIndirect;
186 
187  if (isResonance(id) && (particle->status() == 3 || particle->status() == 22) ){
188  return kDirect;
189  }
190 
191 
192  if (!isIgnored(id) && (isParton(id)))
193  return kNo;
194 
195 
196 
197  unsigned int nMo=particle->numberOfMothers();
198  if (!nMo)
199  return kNo;
200 
201 
202  for(unsigned int i=0;i<nMo;++i){
203  ResonanceState result = fromResonance(invalid,p,dynamic_cast<const reco::GenParticle*>(particle->mother(i)));
204  switch(result) {
205  case kNo:
206  break;
207  case kDirect:
208  if (dynamic_cast<const reco::GenParticle*>(particle->mother(i))->pdgId()==id || isResonance(id))
209  return kDirect;
210  if(!isExcludedFromResonance(id))
211  break;
212  case kIndirect:
213  return kIndirect;
214  }
215  }
216 return kNo;
217 }
218 
219 
221  const ParticleVector &p,
222  const reco::GenParticle *particle) const {
223  return testPartonChildren(invalid, p, particle) > 0;
224 }
225 
226 //######################################################
227 //function NEEDED and called per EVENT by FRAMEWORK:
229 
230 
231  std::auto_ptr<reco::GenParticleRefVector> selected_ (new reco::GenParticleRefVector);
232 
234  evt.getByToken(input_genpartcoll_token_, genParticles );
235 
236  std::map<const reco::GenParticle*,size_t> particlePtrIdxMap;
237  ParticleVector particles;
238  for (reco::GenParticleCollection::const_iterator iter=genParticles->begin();iter!=genParticles->end();++iter){
239  particles.push_back(&*iter);
240  particlePtrIdxMap[&*iter] = (iter - genParticles->begin());
241  }
242 
243  std::sort(particles.begin(), particles.end());
244  unsigned int size = particles.size();
245 
246  ParticleBitmap selected(size, false);
247  ParticleBitmap invalid(size, false);
248 
249  for(unsigned int i = 0; i < size; i++) {
250  const reco::GenParticle *particle = particles[i];
251  if (invalid[i])
252  continue;
253  if (particle->status() == 1)
254  selected[i] = true;
255  if (partonicFinalState && isParton(particle->pdgId())) {
256 
257  if (particle->numberOfDaughters()==0 &&
258  particle->status() != 1) {
259  // some brokenness in event...
260  invalid[i] = true;
261  }
262  else if (!hasPartonChildren(invalid, particles,
263  particle)) {
264  selected[i] = true;
265  invalidateTree(invalid, particles,particle); //this?!?
266  }
267  }
268 
269  }
270 
271  for(size_t idx = 0; idx < size; ++idx){
272  const reco::GenParticle *particle = particles[idx];
273  if (!selected[idx] || invalid[idx]){
274  continue;
275  }
276 
277  if (excludeResonances &&
278  fromResonance(invalid, particles, particle)) {
279  invalid[idx] = true;
280  //cout<<"[INPUTSELECTOR] Invalidates FROM RESONANCE!: ["<<setw(4)<<idx<<"] "<<particle->pdgId()<<" "<<particle->pt()<<endl;
281  continue;
282  }
283 
284  if (isIgnored(particle->pdgId())){
285  continue;
286  }
287 
288 
289  if (particle->pt() >= ptMin){
290  edm::Ref<reco::GenParticleCollection> particleRef(genParticles,particlePtrIdxMap[particle]);
291  selected_->push_back(particleRef);
292  //cout<<"Finally we have: ["<<setw(4)<<idx<<"] "<<setw(4)<<particle->pdgId()<<" "<<particle->pt()<<endl;
293  }
294  }
295  evt.put(selected_);
296 }
297 
298 
299 
300 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool hasPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const reco::GenParticle *particle) const
std::vector< const reco::GenParticle * > ParticleVector
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual int pdgId() const GCC11_FINAL
PDG identifier.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int testPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const reco::GenParticle *particle) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
double npart
Definition: HydjetWrapper.h:44
std::vector< unsigned int > excludeFromResonancePids
virtual void produce(edm::Event &evt, const edm::EventSetup &evtSetup)
void setIgnoredParticles(const std::vector< unsigned int > &particleIDs)
edm::EDGetTokenT< reco::GenParticleCollection > input_genpartcoll_token_
virtual int status() const GCC11_FINAL
status word
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
static void invalidateTree(InputGenJetsParticleSelector::ParticleBitmap &invalid, const InputGenJetsParticleSelector::ParticleVector &p, const reco::GenParticle *particle)
virtual size_t numberOfMothers() const
number of mothers
virtual size_t numberOfDaughters() const
number of daughters
tuple result
Definition: query.py:137
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p, const reco::GenParticle *particle)
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...
std::vector< unsigned int > ignoreParticleIDs
void setExcludeFromResonancePids(const std::vector< unsigned int > &particleIDs)
virtual float pt() const GCC11_FINAL
transverse momentum
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
tuple size
Write out results.
ResonanceState fromResonance(ParticleBitmap &invalid, const ParticleVector &p, const reco::GenParticle *particle) const