CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TopDecaySubset.cc
Go to the documentation of this file.
3 
6 
9  srcToken_ (consumes<reco::GenParticleCollection>(cfg.getParameter<edm::InputTag>("src"))),
10  addRadiation_(cfg.getParameter<bool>("addRadiation")),
11  showerModel_ (kStart)
12 {
13  // mapping of the corresponding fillMode; see FillMode
14  // enumerator of TopDecaySubset for available modes
15  std::string mode = cfg.getParameter<std::string>( "fillMode" );
16  if(mode=="kME")
17  fillMode_=kME;
18  else if(mode=="kStable")
20  else
21  throw cms::Exception("Configuration") << mode << " is not a supported FillMode!\n";
22 
23  // produces a set of GenParticles following
24  // the decay branch of top quarks to the first level of
25  // stable decay products
26  produces<reco::GenParticleCollection>();
27 }
28 
31 {
32 }
33 
35 void
37 {
38  // create target vector
39  std::auto_ptr<reco::GenParticleCollection> target( new reco::GenParticleCollection );
40 
41  // get source collection
43  event.getByToken(srcToken_, src);
44 
45  // find top quarks in list of input particles
46  std::vector<const reco::GenParticle*> tops = findTops(*src);
47 
48  // determine shower model (only in first event)
49  if(showerModel_==kStart)
51 
52  if(showerModel_!=kNone) {
53  // check sanity of W bosons
54  checkWBosons(tops);
55 
56  // get ref product from the event
57  const reco::GenParticleRefProd ref = event.getRefBeforePut<reco::GenParticleCollection>();
58  // clear existing refs
60  // fill output
61  fillListing(tops, *target);
62  // fill references
63  fillReferences(ref, *target);
64  }
65 
66  // write vectors to the event
67  event.put(target);
68 }
69 
71 std::vector<const reco::GenParticle*>
73 {
74  std::vector<const reco::GenParticle*> tops;
75  for(reco::GenParticleCollection::const_iterator t=parts.begin(); t!=parts.end(); ++t){
76  if( std::abs(t->pdgId())==TopDecayID::tID && t->status()==TopDecayID::unfrag )
77  tops.push_back( &(*t) );
78  }
79  return tops;
80 }
81 
84 TopDecaySubset::checkShowerModel(const std::vector<const reco::GenParticle*>& tops) const
85 {
86  for(std::vector<const reco::GenParticle*>::const_iterator it=tops.begin(); it!=tops.end(); ++it){
87  const reco::GenParticle* top = *it;
88  // check for kHerwig type showers: here the status 3 top quark will
89  // have a single status 2 top quark as daughter, which has again 3
90  // or more status 2 daughters
91  if( top->numberOfDaughters()==1){
92  if( top->begin()->pdgId()==top->pdgId() && top->begin()->status()==TopDecayID::stable && top->begin()->numberOfDaughters()>1)
93  return kHerwig;
94  }
95  // check for kPythia type showers: here the status 3 top quark will
96  // have all decay products and a status 2 top quark as daughters
97  // the status 2 top quark will be w/o further daughters
98  if( top->numberOfDaughters()>1 ){
99  bool containsWBoson=false, containsQuarkDaughter=false;
100  for(reco::GenParticle::const_iterator td=top->begin(); td!=top->end(); ++td){
101  if( std::abs(td->pdgId ()) <TopDecayID::tID ) containsQuarkDaughter=true;
102  if( std::abs(td->pdgId ())==TopDecayID::WID ) containsWBoson=true;
103  }
104  if(containsQuarkDaughter && containsWBoson)
105  return kPythia;
106  }
107  }
108  // if neither Herwig nor Pythia like
109  if(tops.size()==0)
110  edm::LogInfo("decayChain")
111  << " Failed to find top quarks in decay chain. Will assume that this a \n"
112  << " non-top sample and produce an empty decaySubset.\n";
113  else
115  " Can not find back any of the supported hadronization models. Models \n"
116  " which are supported are: \n"
117  " Pythia LO(+PS): Top(status 3) --> WBoson(status 3), Quark(status 3)\n"
118  " Herwig NLO(+PS): Top(status 2) --> Top(status 3) --> Top(status 2) \n");
119  return kNone;
120 }
121 
123 void
124 TopDecaySubset::checkWBosons(std::vector<const reco::GenParticle*>& tops) const
125 {
126  unsigned nTops = tops.size();
127  for(std::vector<const reco::GenParticle*>::iterator it=tops.begin(); it!=tops.end();) {
128  const reco::GenParticle* top = *it;
129  bool isContained=false;
130  bool expectedStatus=false;
131  if(showerModel_!=kPythia && top->begin()==top->end())
133  "showerModel_!=kPythia && top->begin()==top->end()\n");
135  td!=((showerModel_==kPythia) ? top->end() : top->begin()->end());
136  ++td){
137  if(std::abs(td->pdgId())==TopDecayID::WID) {
138  isContained=true;
139  if( ((showerModel_==kPythia) ? td->status()==TopDecayID::unfrag : td->status()==TopDecayID::stable) ) {
140  expectedStatus=true;
141  break;
142  }
143  }
144  }
145  if(!expectedStatus) {
146  it=tops.erase(it);
147  if(isContained)
148  edm::LogInfo("decayChain")
149  << " W boson does not have the expected status. This happens, e.g., \n"
150  << " with MC@NLO in the case of additional ttbar pairs from radiated \n"
151  << " gluons. We hope everything is fine, remove the correspdonding top \n"
152  << " quark from our list since it is not part of the primary ttbar pair \n"
153  << " and try to continue. \n";
154  }
155  else
156  it++;
157  }
158  if(tops.size()==0 && nTops!=0)
160  " Did not find a W boson with appropriate status for any of the top \n"
161  " quarks in this event. This means that the hadronization of the W \n"
162  " boson might be screwed up or there is another problem with the \n"
163  " particle listing in this MC sample. Please contact an expert. \n");
164 }
165 
167 void
168 TopDecaySubset::fillListing(const std::vector<const reco::GenParticle*>& tops, reco::GenParticleCollection& target)
169 {
170  unsigned int statusFlag;
171  // determine status flag of the new
172  // particle depending on the FillMode
173  fillMode_ == kME ? statusFlag=3 : statusFlag=2;
174 
175  for(std::vector<const reco::GenParticle*>::const_iterator it=tops.begin(); it!=tops.end(); ++it){
176  const reco::GenParticle* t = *it;
177  // if particle is top or anti-top
178  std::auto_ptr<reco::GenParticle> topPtr( new reco::GenParticle( t->threeCharge(), p4(it, statusFlag), t->vertex(), t->pdgId(), statusFlag, false ) );
179  target.push_back( *topPtr );
180  ++motherPartIdx_;
181  // keep the top index for the map to manage the daughter refs
182  int iTop=motherPartIdx_;
183  std::vector<int> topDaughters;
184  // define the W boson index (to be set later) for the map to
185  // manage the daughter refs
186  int iW = 0;
187  std::vector<int> wDaughters;
188  // sanity check
189  if(showerModel_!=kPythia && t->begin()==t->end())
191  "showerModel_!=kPythia && t->begin()==t->end()\n");
192  //iterate over top daughters
194  if( td->status()==TopDecayID::unfrag && std::abs( td->pdgId() )<=TopDecayID::bID ){
195  // if particle is beauty or other quark
196  std::auto_ptr<reco::GenParticle> bPtr( new reco::GenParticle( td->threeCharge(), p4( td, statusFlag ), td->vertex(), td->pdgId(), statusFlag, false ) );
197  target.push_back( *bPtr );
198  // increment & push index of the top daughter
199  topDaughters.push_back( ++motherPartIdx_ );
200  if(addRadiation_){
202  }
203  }
204  // sanity check
205  if(showerModel_!=kPythia && td->begin()==td->end())
207  "showerModel_!=kPythia && td->begin()==td->end()\n");
209  if( buffer->status()==TopDecayID::unfrag && std::abs( buffer->pdgId() )==TopDecayID::WID ){
210  // if particle is a W boson
211  std::auto_ptr<reco::GenParticle> wPtr( new reco::GenParticle( buffer->threeCharge(), p4( buffer, statusFlag), buffer->vertex(), buffer->pdgId(), statusFlag, true ) );
212  target.push_back( *wPtr );
213  // increment & push index of the top daughter
214  topDaughters.push_back( ++motherPartIdx_ );
215  // keep the W idx for the map
216  iW=motherPartIdx_;
217  if(addRadiation_){
218  addRadiation(motherPartIdx_,buffer,target);
219  }
220  if(showerModel_!=kPythia && buffer->begin()==buffer->end())
222  "showerModel_!=kPythia && buffer->begin()==buffer->end()\n");
223  // iterate over W daughters
224  for(reco::GenParticle::const_iterator wd=((showerModel_==kPythia)?buffer->begin():buffer->begin()->begin()); wd!=((showerModel_==kPythia)?buffer->end():buffer->begin()->end()); ++wd){
225  // make sure the W daughter is of status unfrag and not the W itself
226  if( wd->status()==TopDecayID::unfrag && !(std::abs(wd->pdgId())==TopDecayID::WID) ) {
227  std::auto_ptr<reco::GenParticle> qPtr( new reco::GenParticle( wd->threeCharge(), p4( wd, statusFlag ), wd->vertex(), wd->pdgId(), statusFlag, false) );
228  target.push_back( *qPtr );
229  // increment & push index of the top daughter
230  wDaughters.push_back( ++motherPartIdx_ );
231  if( wd->status()==TopDecayID::unfrag && std::abs( wd->pdgId() )==TopDecayID::tauID ){
232  // add tau daughters if the particle is a tau pass
233  // the daughter of the tau which is of status 2
234  //addDaughters(motherPartIdx_, wd->begin(), target);
235  // add tau daughters if the particle is a tau pass
236  // the tau itself, which may add a tau daughter of
237  // of status 2 to the listing
238  addDaughters(motherPartIdx_,wd,target);
239  }
240  }
241  }
242  }
243  if(addRadiation_ && buffer->status()==TopDecayID::stable && ( buffer->pdgId()==TopDecayID::glueID || std::abs(buffer->pdgId())<TopDecayID::bID)){
244  // collect additional radiation from the top
245  std::auto_ptr<reco::GenParticle> radPtr( new reco::GenParticle( buffer->threeCharge(), buffer->p4(), buffer->vertex(), buffer->pdgId(), statusFlag, false ) );
246  target.push_back( *radPtr );
247  // increment & push index of the top daughter
248  topDaughters.push_back( ++motherPartIdx_ );
249  }
250  }
251  // add potential sisters of the top quark;
252  // only for top to prevent double counting
253  if(t->numberOfMothers()>0 && t->pdgId()==TopDecayID::tID){
254  for(reco::GenParticle::const_iterator ts = t->mother()->begin(); ts!=t->mother()->end(); ++ts){
255  // loop over all daughters of the top mother i.e.
256  // the two top quarks and their potential sisters
257  if(std::abs(ts->pdgId())!=t->pdgId() && ts->pdgId()!=t->mother()->pdgId()){
258  // add all further particles but the two top quarks and potential
259  // cases where the mother of the top has itself as daughter
260  reco::GenParticle* cand = new reco::GenParticle( ts->threeCharge(), ts->p4(), ts->vertex(), ts->pdgId(), ts->status(), false );
261  std::auto_ptr<reco::GenParticle> sPtr( cand );
262  target.push_back( *sPtr );
263  if(ts->begin()!=ts->end()){
264  // in case the sister has daughters increment
265  // and add the first generation of daughters
266  addDaughters(++motherPartIdx_,ts->begin(),target,false);
267  }
268  }
269  }
270  }
271  // fill the map for the administration
272  // of daughter indices
273  refs_[ iTop ] = topDaughters;
274  refs_[ iW ] = wDaughters;
275  }
276 }
277 
280 TopDecaySubset::p4(const std::vector<const reco::GenParticle*>::const_iterator top, int statusFlag)
281 {
282  // calculate the four vector for top/anti-top quarks from
283  // the W boson and the b quark plain or including all
284  // additional radiation depending on switch 'plain'
285  if(statusFlag==TopDecayID::unfrag){
286  // return 4 momentum as it is
287  return (*top)->p4();
288  }
290  for(reco::GenParticle::const_iterator p=(*top)->begin(); p!=(*top)->end(); ++p){
291  if( p->status() == TopDecayID::unfrag ){
292  // descend by one level for each
293  // status 3 particle on the way
294  vec+=p4( p, statusFlag );
295  }
296  else{
297  if( std::abs((*top)->pdgId())==TopDecayID::WID ){
298  // in case of a W boson skip the status
299  // 2 particle to prevent double counting
300  if( std::abs(p->pdgId())!=TopDecayID::WID )
301  vec+=p->p4();
302  }
303  else{
304  // add all four vectors for each stable
305  // particle (status 1 or 2) on the way
306  vec+=p->p4();
307  if( vec.mass()-(*top)->mass()>0 ){
308  // continue adding up gluons and qqbar pairs on the top
309  // line untill the nominal top mass is reached; then
310  // break in order to prevent picking up virtualities
311  break;
312  }
313  }
314  }
315  }
316  return vec;
317 }
318 
322 {
323  // calculate the four vector for all top daughters from
324  // their daughters including additional radiation
325  if(statusFlag==TopDecayID::unfrag){
326  // return 4 momentum as it is
327  return part->p4();
328  }
330  for(reco::GenParticle::const_iterator p=part->begin(); p!=part->end(); ++p){
331  if( p->status()<=TopDecayID::stable && std::abs(p->pdgId ())==TopDecayID::WID){
332  vec=p->p4();
333  }
334  else{
335  if(p->status()<=TopDecayID::stable){
336  // sum up the p4 of all stable particles
337  // (of status 1 or 2)
338  vec+=p->p4();
339  }
340  else{
341  if( p->status()==TopDecayID::unfrag){
342  // if the particle is unfragmented (i.e.
343  // status 3) descend by one level
344  vec+=p4(p, statusFlag);
345  }
346  }
347  }
348  }
349  return vec;
350 }
351 
353 void
355 {
356  std::vector<int> daughters;
357  int idxBuffer = idx;
358  for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
359  if(daughter->status()<=TopDecayID::stable && daughter->pdgId()!=part->pdgId()){
360  // skip comment lines and make sure that
361  // the particle is not double counted as
362  // dasughter of itself
363  std::auto_ptr<reco::GenParticle> ptr( new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
364  target.push_back( *ptr );
365  daughters.push_back( ++idx ); //push index of daughter
366  }
367  }
368  if(daughters.size()) {
369  refs_[ idxBuffer ] = daughters;
370  }
371 }
372 
374 void
376 {
377  std::vector<int> daughters;
378  int idxBuffer = idx;
379  for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
380  std::auto_ptr<reco::GenParticle> ptr( new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
381  target.push_back( *ptr );
382  // increment & push index of daughter
383  daughters.push_back( ++idx );
384  // continue recursively if desired
385  if(recursive){
386  addDaughters(idx,daughter,target);
387  }
388  }
389  if(daughters.size()) {
390  refs_[ idxBuffer ] = daughters;
391  }
392 }
393 
395 void
397 {
398  // clear vector of references
399  refs_.clear();
400  // set idx for mother particles to a start value
401  // of -1 (the first entry will raise it to 0)
402  motherPartIdx_=-1;
403 }
404 
406 void
408 {
409  int idx=0;
410  for(reco::GenParticleCollection::iterator p=sel.begin(); p!=sel.end(); ++p, ++idx){
411  //find daughter reference vectors in refs_ and add daughters
412  std::map<int, std::vector<int> >::const_iterator daughters=refs_.find( idx );
413  if( daughters!=refs_.end() ){
414  for(std::vector<int>::const_iterator daughter = daughters->second.begin();
415  daughter!=daughters->second.end(); ++daughter){
416  reco::GenParticle* part = dynamic_cast<reco::GenParticle* > (&(*p));
417  if(part == 0){
418  throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
419  }
420  part->addDaughter( reco::GenParticleRef(ref, *daughter) );
421  sel[*daughter].addMother( reco::GenParticleRef(ref, idx) );
422  }
423  }
424  }
425 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
virtual int pdgId() const
PDG identifier.
static bool isContained(const std::vector< unsigned int > &list, int id)
Definition: JetInput.cc:89
tuple cfg
Definition: looper.py:237
void addDaughters(int &idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection &target, bool recursive=true)
recursively fill vector for all further decay particles of a given particle
static const int bID
Definition: TopGenEvent.h:14
ShowerModel checkShowerModel(const std::vector< const reco::GenParticle * > &tops) const
check the decay chain for the used shower model
ShowerModel showerModel_
parton shower mode (filled in checkShowerModel)
virtual void produce(edm::Event &event, const edm::EventSetup &setup)
write output into the event
virtual const Point & vertex() const
vertex position (overwritten by PF...)
static const int glueID
Definition: TopGenEvent.h:15
void addRadiation(int &idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection &target)
fill vector including all radiations from quarks originating from W/top
static const int unfrag
Definition: TopGenEvent.h:12
virtual int threeCharge() const =0
electric charge
virtual int status() const =0
status word
static const int stable
Definition: TopGenEvent.h:11
std::vector< const reco::GenParticle * > findTops(const reco::GenParticleCollection &parts)
find top quarks in list of input particles
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
void clearReferences()
clear references
static const int tID
Definition: TopGenEvent.h:13
static const int tauID
Definition: TopGenEvent.h:21
virtual size_type numberOfDaughters() const =0
number of daughters
virtual const_iterator end() const =0
last daughter const_iterator
virtual const_iterator begin() const
first daughter const_iterator
virtual size_t numberOfMothers() const
number of mothers
virtual size_t numberOfDaughters() const
number of daughters
reco::Particle::LorentzVector p4(const std::vector< const reco::GenParticle * >::const_iterator top, int statusFlag)
calculate lorentz vector from input (dedicated to top reconstruction)
void checkWBosons(std::vector< const reco::GenParticle * > &tops) const
check whether W bosons are contained in the original gen particle listing
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool addRadiation_
add radiation or not?
ShowerModel
classification of potential shower types
virtual const Point & vertex() const =0
vertex position
virtual int threeCharge() const
electric charge
edm::EDGetTokenT< reco::GenParticleCollection > srcToken_
input tag for the genParticle source
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
void fillReferences(const reco::GenParticleRefProd &refProd, reco::GenParticleCollection &target)
fill references for output vector
virtual int pdgId() const =0
PDG identifier.
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
part
Definition: HCALResponse.h:20
TopDecaySubset(const edm::ParameterSet &cfg)
default constructor
static const int WID
Definition: TopGenEvent.h:18
std::map< int, std::vector< int > > refs_
management of daughter indices for fillRefs
void fillListing(const std::vector< const reco::GenParticle * > &tops, reco::GenParticleCollection &target)
fill output vector for full decay chain
FillMode fillMode_
virtual const_iterator begin() const =0
first daughter const_iterator
volatile std::atomic< bool > shutdown_flag false
~TopDecaySubset()
default destructor
virtual const_iterator end() const
last daughter const_iterator
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector