CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/TopQuarkAnalysis/TopEventProducers/src/TopDecaySubset.cc

Go to the documentation of this file.
00001 #include "FWCore/Utilities/interface/EDMException.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 
00004 #include "AnalysisDataFormats/TopObjects/interface/TopGenEvent.h"
00005 #include "TopQuarkAnalysis/TopEventProducers/interface/TopDecaySubset.h"
00006 
00008 TopDecaySubset::TopDecaySubset(const edm::ParameterSet& cfg): 
00009   src_         (cfg.getParameter<edm::InputTag>("src")),
00010   addRadiation_(cfg.getParameter<bool>("addRadiation")),
00011   showerModel_ (kStart)
00012 {
00013   // mapping of the corresponding fillMode; see FillMode 
00014   // enumerator of TopDecaySubset for available modes
00015   std::string mode = cfg.getParameter<std::string>( "fillMode" );
00016   if(mode=="kME")
00017     fillMode_=kME;
00018   else if(mode=="kStable")
00019     fillMode_=kStable;
00020   else
00021     throw cms::Exception("Configuration") << mode << " is not a supported FillMode!\n";
00022 
00023   // produces a set of GenParticles following
00024   // the decay branch of top quarks to the first level of 
00025   // stable decay products
00026   produces<reco::GenParticleCollection>(); 
00027 }
00028 
00030 TopDecaySubset::~TopDecaySubset()
00031 {
00032 }
00033 
00035 void
00036 TopDecaySubset::produce(edm::Event& event, const edm::EventSetup& setup)
00037 {     
00038   // get source collection
00039   edm::Handle<reco::GenParticleCollection> src;
00040   event.getByLabel(src_, src);
00041 
00042   // find top quarks in list of input particles
00043   std::vector<const reco::GenParticle*> tops = findTops(*src);
00044 
00045   // determine shower model (only in first event)
00046   if(showerModel_==kStart)
00047     showerModel_=checkShowerModel(tops);
00048 
00049   // check sanity of W bosons
00050   checkWBosons(tops);
00051 
00052   // create target vector
00053   std::auto_ptr<reco::GenParticleCollection> target( new reco::GenParticleCollection );
00054   // get ref product from the event
00055   const reco::GenParticleRefProd ref = event.getRefBeforePut<reco::GenParticleCollection>(); 
00056   // clear existing refs
00057   clearReferences();  
00058   // fill output
00059   fillListing(tops, *target);
00060   // fill references
00061   fillReferences(ref, *target);
00062 
00063   // write vectors to the event
00064   event.put(target);
00065 }
00066 
00068 std::vector<const reco::GenParticle*>
00069 TopDecaySubset::findTops(const reco::GenParticleCollection& parts)
00070 {
00071   std::vector<const reco::GenParticle*> tops;
00072   for(reco::GenParticleCollection::const_iterator t=parts.begin(); t!=parts.end(); ++t){
00073     if( std::abs(t->pdgId())==TopDecayID::tID && t->status()==TopDecayID::unfrag )
00074       tops.push_back( &(*t) );
00075   }
00076   return tops;
00077 }
00078 
00080 TopDecaySubset::ShowerModel
00081 TopDecaySubset::checkShowerModel(const std::vector<const reco::GenParticle*>& tops) const
00082 {
00083   for(std::vector<const reco::GenParticle*>::const_iterator it=tops.begin(); it!=tops.end(); ++it){
00084     const reco::GenParticle* top = *it;
00085     // check for kHerwig type showers: here the status 3 top quark will 
00086     // have a single status 2 top quark as daughter, which has again 3 
00087     // or more status 2 daughters
00088     if( top->numberOfDaughters()==1){
00089       if( top->begin()->pdgId()==top->pdgId() && top->begin()->status()==TopDecayID::stable && top->begin()->numberOfDaughters()>1)
00090         return kHerwig;
00091     }
00092     // check for kPythia type showers: here the status 3 top quark will 
00093     // have all decay products and a status 2 top quark as daughters
00094     // the status 2 top quark will be w/o further daughters
00095     if( top->numberOfDaughters()>1 ){
00096       bool containsWBoson=false, containsQuarkDaughter=false;
00097       for(reco::GenParticle::const_iterator td=top->begin(); td!=top->end(); ++td){
00098         if( std::abs(td->pdgId ()) <TopDecayID::tID ) containsQuarkDaughter=true;
00099         if( std::abs(td->pdgId ())==TopDecayID::WID ) containsWBoson=true;
00100       }
00101       if(containsQuarkDaughter && containsWBoson)
00102         return kPythia;
00103     }
00104   }
00105   // if neither Herwig nor Pythia like
00106   if(tops.size()==0)
00107     edm::LogWarning("decayChain")
00108       << " Failed to find top quarks in decay chain. Will assume that this a \n"
00109       << " non-top sample and produce an empty decaySubset.\n";
00110   else
00111     throw edm::Exception(edm::errors::LogicError,
00112                          " Can not find back any of the supported hadronization models. Models \n"
00113                          " which are supported are:                                            \n"
00114                          " Pythia  LO(+PS): Top(status 3) --> WBoson(status 3), Quark(status 3)\n"
00115                          " Herwig NLO(+PS): Top(status 2) --> Top(status 3) --> Top(status 2)  \n");
00116   return kNone;
00117 }
00118 
00120 void 
00121 TopDecaySubset::checkWBosons(std::vector<const reco::GenParticle*>& tops) const
00122 {
00123   unsigned nTops = tops.size();
00124   for(std::vector<const reco::GenParticle*>::iterator it=tops.begin(); it!=tops.end();) {
00125     const reco::GenParticle* top = *it;
00126     bool isContained=false;
00127     bool expectedStatus=false;
00128     for(reco::GenParticle::const_iterator td=((showerModel_==kPythia) ? top->begin() : top->begin()->begin());
00129         td!=((showerModel_==kPythia) ? top->end() : top->begin()->end());
00130         ++td){
00131       if(std::abs(td->pdgId())==TopDecayID::WID) {
00132         isContained=true;
00133         if( ((showerModel_==kPythia) ? td->status()==TopDecayID::unfrag : td->status()==TopDecayID::stable) ) { 
00134           expectedStatus=true;
00135           break;
00136         }
00137       }
00138     }
00139     if(!expectedStatus) {
00140       it=tops.erase(it);
00141       if(isContained)
00142         edm::LogWarning("decayChain")
00143           << " W boson does not have the expected status. This happens, e.g.,      \n"
00144           << " with MC@NLO in the case of additional ttbar pairs from radiated     \n"
00145           << " gluons. We hope everything is fine, remove the correspdonding top   \n"
00146           << " quark from our list since it is not part of the primary ttbar pair  \n"
00147           << " and try to continue.                                                \n";      
00148     }
00149     else
00150       it++;
00151   }
00152   if(tops.size()==0 && nTops!=0)
00153     throw edm::Exception(edm::errors::LogicError,
00154                          " Did not find a W boson with appropriate status for any of the top   \n"
00155                          " quarks in this event. This means that the hadronization of the W    \n"
00156                          " boson might be screwed up or there is another problem with the      \n"
00157                          " particle listing in this MC sample. Please contact an expert.       \n");
00158 }
00159 
00161 void 
00162 TopDecaySubset::fillListing(const std::vector<const reco::GenParticle*>& tops, reco::GenParticleCollection& target)
00163 {
00164   unsigned int statusFlag;
00165   // determine status flag of the new 
00166   // particle depending on the FillMode
00167   fillMode_ == kME ? statusFlag=3 : statusFlag=2;
00168 
00169   for(std::vector<const reco::GenParticle*>::const_iterator it=tops.begin(); it!=tops.end(); ++it){
00170     const reco::GenParticle* t = *it;
00171     // if particle is top or anti-top 
00172     std::auto_ptr<reco::GenParticle> topPtr( new reco::GenParticle( t->threeCharge(), p4(it, statusFlag), t->vertex(), t->pdgId(), statusFlag, false ) );
00173     target.push_back( *topPtr );
00174     ++motherPartIdx_;
00175     // keep the top index for the map to manage the daughter refs
00176     int iTop=motherPartIdx_; 
00177     std::vector<int> topDaughters;
00178     // define the W boson index (to be set later) for the map to 
00179     // manage the daughter refs
00180     int iW = 0;
00181     std::vector<int> wDaughters;
00182     //iterate over top daughters
00183     for(reco::GenParticle::const_iterator td=((showerModel_==kPythia)?t->begin():t->begin()->begin()); td!=((showerModel_==kPythia)?t->end():t->begin()->end()); ++td){
00184       if( td->status()==TopDecayID::unfrag && std::abs( td->pdgId() )<=TopDecayID::bID ){ 
00185         // if particle is beauty or other quark 
00186         std::auto_ptr<reco::GenParticle> bPtr( new reco::GenParticle( td->threeCharge(), p4( td, statusFlag ), td->vertex(), td->pdgId(), statusFlag, false ) );
00187         target.push_back( *bPtr );        
00188         // increment & push index of the top daughter
00189         topDaughters.push_back( ++motherPartIdx_ ); 
00190         if(addRadiation_){
00191           addRadiation(motherPartIdx_,td,target); 
00192         }
00193       }
00194       reco::GenParticle::const_iterator buffer = (showerModel_==kPythia)?td:td->begin();
00195       if( buffer->status()==TopDecayID::unfrag && std::abs( buffer->pdgId() )==TopDecayID::WID ){ 
00196         // if particle is a W boson
00197         std::auto_ptr<reco::GenParticle> wPtr(  new reco::GenParticle( buffer->threeCharge(), p4( buffer, statusFlag), buffer->vertex(), buffer->pdgId(), statusFlag, true ) );
00198         target.push_back( *wPtr );
00199         // increment & push index of the top daughter
00200         topDaughters.push_back( ++motherPartIdx_ );
00201         // keep the W idx for the map
00202         iW=motherPartIdx_; 
00203         if(addRadiation_){
00204           addRadiation(motherPartIdx_,buffer,target); 
00205         }
00206         // iterate over W daughters
00207         for(reco::GenParticle::const_iterator wd=((showerModel_==kPythia)?buffer->begin():buffer->begin()->begin()); wd!=((showerModel_==kPythia)?buffer->end():buffer->begin()->end()); ++wd){
00208           // make sure the W daughter is of status unfrag and not the W itself
00209           if( wd->status()==TopDecayID::unfrag && !(std::abs(wd->pdgId())==TopDecayID::WID) ) {
00210             std::auto_ptr<reco::GenParticle> qPtr( new reco::GenParticle( wd->threeCharge(), p4( wd, statusFlag ), wd->vertex(), wd->pdgId(), statusFlag, false) );
00211             target.push_back( *qPtr );
00212             // increment & push index of the top daughter
00213             wDaughters.push_back( ++motherPartIdx_ );
00214             if( wd->status()==TopDecayID::unfrag && std::abs( wd->pdgId() )==TopDecayID::tauID ){ 
00215               // add tau daughters if the particle is a tau pass
00216               // the daughter of the tau which is of status 2
00217               //addDaughters(motherPartIdx_, wd->begin(), target); 
00218               // add tau daughters if the particle is a tau pass
00219               // the tau itself, which may add a tau daughter of
00220               // of status 2 to the listing
00221               addDaughters(motherPartIdx_,wd,target); 
00222             }
00223           } 
00224         }
00225       }
00226       if(addRadiation_ && buffer->status()==TopDecayID::stable && ( buffer->pdgId()==TopDecayID::glueID || std::abs(buffer->pdgId())<TopDecayID::bID)){
00227         // collect additional radiation from the top 
00228         std::auto_ptr<reco::GenParticle> radPtr( new reco::GenParticle( buffer->threeCharge(), buffer->p4(), buffer->vertex(), buffer->pdgId(), statusFlag, false ) );
00229         target.push_back( *radPtr );    
00230         // increment & push index of the top daughter  
00231         topDaughters.push_back( ++motherPartIdx_ );
00232       }
00233     }
00234     // add potential sisters of the top quark;
00235     // only for top to prevent double counting
00236     if(t->numberOfMothers()>0 && t->pdgId()==TopDecayID::tID){
00237       for(reco::GenParticle::const_iterator ts = t->mother()->begin(); ts!=t->mother()->end(); ++ts){
00238         // loop over all daughters of the top mother i.e.
00239         // the two top quarks and their potential sisters
00240         if(std::abs(ts->pdgId())!=t->pdgId() && ts->pdgId()!=t->mother()->pdgId()){
00241           // add all further particles but the two top quarks and potential 
00242           // cases where the mother of the top has itself as daughter
00243           reco::GenParticle* cand = new reco::GenParticle( ts->threeCharge(), ts->p4(), ts->vertex(), ts->pdgId(), ts->status(), false );
00244           std::auto_ptr<reco::GenParticle> sPtr( cand );
00245           target.push_back( *sPtr );
00246           if(ts->begin()!=ts->end()){ 
00247             // in case the sister has daughters increment
00248             // and add the first generation of daughters
00249             addDaughters(++motherPartIdx_,ts->begin(),target,false);
00250           }
00251         }
00252       }
00253     }
00254     // fill the map for the administration 
00255     // of daughter indices
00256     refs_[ iTop ] = topDaughters;
00257     refs_[ iW   ] = wDaughters; 
00258   }
00259 }
00260 
00262 reco::Particle::LorentzVector 
00263 TopDecaySubset::p4(const std::vector<const reco::GenParticle*>::const_iterator top, int statusFlag)
00264 {
00265   // calculate the four vector for top/anti-top quarks from 
00266   // the W boson and the b quark plain or including all 
00267   // additional radiation depending on switch 'plain'
00268   if(statusFlag==TopDecayID::unfrag){
00269     // return 4 momentum as it is
00270     return (*top)->p4();
00271   }
00272   reco::Particle::LorentzVector vec;
00273   for(reco::GenParticle::const_iterator p=(*top)->begin(); p!=(*top)->end(); ++p){
00274     if( p->status() == TopDecayID::unfrag ){
00275       // descend by one level for each
00276       // status 3 particle on the way
00277       vec+=p4( p, statusFlag );
00278     }
00279     else{ 
00280       if( std::abs((*top)->pdgId())==TopDecayID::WID ){
00281         // in case of a W boson skip the status 
00282         // 2 particle to prevent double counting
00283         if( std::abs(p->pdgId())!=TopDecayID::WID ) 
00284           vec+=p->p4();
00285       } 
00286       else{
00287         // add all four vectors for each stable
00288         // particle (status 1 or 2) on the way 
00289         vec+=p->p4();
00290         if( vec.mass()-(*top)->mass()>0 ){
00291           // continue adding up gluons and qqbar pairs on the top 
00292           // line untill the nominal top mass is reached; then 
00293           // break in order to prevent picking up virtualities
00294           break;
00295         }
00296       }
00297     }
00298   }
00299   return vec;
00300 }
00301 
00303 reco::Particle::LorentzVector 
00304 TopDecaySubset::p4(const reco::GenParticle::const_iterator part, int statusFlag)
00305 {
00306   // calculate the four vector for all top daughters from 
00307   // their daughters including additional radiation 
00308   if(statusFlag==TopDecayID::unfrag){
00309     // return 4 momentum as it is
00310     return part->p4();
00311   }
00312   reco::Particle::LorentzVector vec;
00313   for(reco::GenParticle::const_iterator p=part->begin(); p!=part->end(); ++p){
00314     if( p->status()<=TopDecayID::stable && std::abs(p->pdgId ())==TopDecayID::WID){
00315       vec=p->p4();
00316     }
00317     else{
00318       if(p->status()<=TopDecayID::stable){
00319         // sum up the p4 of all stable particles 
00320         // (of status 1 or 2)
00321         vec+=p->p4();
00322       }
00323       else{
00324         if( p->status()==TopDecayID::unfrag){
00325           // if the particle is unfragmented (i.e.
00326           // status 3) descend by one level
00327           vec+=p4(p, statusFlag);   
00328         }
00329       }
00330     }
00331   }
00332   return vec;
00333 }
00334 
00336 void 
00337 TopDecaySubset::addRadiation(int& idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection& target)
00338 {
00339   std::vector<int> daughters;
00340   int idxBuffer = idx;
00341   for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
00342     if(daughter->status()<=TopDecayID::stable && daughter->pdgId()!=part->pdgId()){
00343       // skip comment lines and make sure that
00344       // the particle is not double counted as 
00345       // dasughter of itself
00346       std::auto_ptr<reco::GenParticle> ptr(  new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
00347       target.push_back( *ptr );
00348       daughters.push_back( ++idx ); //push index of daughter
00349     }
00350   }  
00351   if(daughters.size()) {
00352     refs_[ idxBuffer ] = daughters;
00353   }
00354 }
00355 
00357 void 
00358 TopDecaySubset::addDaughters(int& idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection& target, bool recursive)
00359 {
00360   std::vector<int> daughters;
00361   int idxBuffer = idx;
00362   for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
00363       std::auto_ptr<reco::GenParticle> ptr( new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
00364       target.push_back( *ptr );
00365       // increment & push index of daughter
00366       daughters.push_back( ++idx );
00367       // continue recursively if desired
00368       if(recursive){
00369         addDaughters(idx,daughter,target);  
00370       }
00371   }  
00372   if(daughters.size()) {
00373     refs_[ idxBuffer ] = daughters;
00374   }
00375 }
00376 
00378 void
00379 TopDecaySubset::clearReferences()
00380 {
00381   // clear vector of references 
00382   refs_.clear();  
00383   // set idx for mother particles to a start value
00384   // of -1 (the first entry will raise it to 0)
00385   motherPartIdx_=-1;
00386 }
00387 
00389 void 
00390 TopDecaySubset::fillReferences(const reco::GenParticleRefProd& ref, reco::GenParticleCollection& sel)
00391 { 
00392   int idx=0;
00393   for(reco::GenParticleCollection::iterator p=sel.begin(); p!=sel.end(); ++p, ++idx){
00394     //find daughter reference vectors in refs_ and add daughters
00395     std::map<int, std::vector<int> >::const_iterator daughters=refs_.find( idx );
00396     if( daughters!=refs_.end() ){
00397       for(std::vector<int>::const_iterator daughter = daughters->second.begin(); 
00398           daughter!=daughters->second.end(); ++daughter){
00399         reco::GenParticle* part = dynamic_cast<reco::GenParticle* > (&(*p));
00400         if(part == 0){
00401          throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
00402         }
00403         part->addDaughter( reco::GenParticleRef(ref, *daughter) );
00404         sel[*daughter].addMother( reco::GenParticleRef(ref, idx) );
00405       }
00406     }
00407   }
00408 }