CMS 3D CMS Logo

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