CMS 3D CMS Logo

TopDecaySubset Class Reference

#include <TopQuarkAnalysis/TopEventProducers/interface/TopDecaySubset.h>

Inheritance diagram for TopDecaySubset:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Types

enum  FillMode { kStable, kME }
 supported modes to fill the new vectors of gen particles More...

Public Member Functions

virtual void produce (edm::Event &, const edm::EventSetup &)
 write TopDecaySubset into the event
 TopDecaySubset (const edm::ParameterSet &)
 default constructor
 ~TopDecaySubset ()
 default destructor

Private Member Functions

void addDaughters (int &idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection &target, bool recursive=true)
 recursively fill vector for all further dacay particles of a given particle
void addRadiation (int &idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection &target)
 fill vector including all radiations from quarks originating from W/top
void clearReferences ()
 clear references
void fillOutput (edm::Event &evt, const reco::GenParticleCollection &src, FillMode mode)
 fill output collection depending on whether the W boson is contained in the generator listing or not
void fillReferences (const reco::GenParticleRefProd &refProd, reco::GenParticleCollection &target)
 fill references for output vector
void fromFullListing (const reco::GenParticleCollection &src, reco::GenParticleCollection &target, const int &partId, FillMode mode)
 fill output vector with full decay chain with intermediate W bosons
void fromTruncListing (const reco::GenParticleCollection &src, reco::GenParticleCollection &target, const int &partId, FillMode mode)
 fill output vector with full decay chain w/o intermediate W bosons
reco::Particle::LorentzVector p4 (const reco::GenParticle::const_iterator part, int statusFlag)
 calculate lorentz vector from input
reco::Particle::LorentzVector p4 (const std::vector< reco::GenParticle >::const_iterator top, int statusFlag)
 calculate lorentz vector from input (dedicated to top reconstruction)
void printSource (const reco::GenParticleCollection &src)
 print the whole listing if particle with pdgId is contained in the top decay chain
void printTarget (reco::GenParticleCollection &target)
 print the whole decay chain if particle with pdgId is contained in the top decay chain
bool wInDecayChain (const reco::GenParticleCollection &src, const int &partId)
 check whether the W boson is contained in the original gen particle listing or not

Private Attributes

bool addRadiatedGluons_
 add radiated gluons or not?
int motherPartIdx_
 index in new evt listing of parts with daughters; has to be set to -1 in produce to deliver consistent results!
std::map< int, std::vector< int > > refs_
 management of daughter indices for fillRefs
edm::InputTag src_
 input tag for the genParticle source


Detailed Description

Definition at line 15 of file TopDecaySubset.h.


Member Enumeration Documentation

enum TopDecaySubset::FillMode

supported modes to fill the new vectors of gen particles

Enumerator:
kStable 
kME 

Definition at line 21 of file TopDecaySubset.h.

00021 {kStable, kME};


Constructor & Destructor Documentation

TopDecaySubset::TopDecaySubset ( const edm::ParameterSet cfg  )  [explicit]

default constructor

Definition at line 14 of file TopDecaySubset.cc.

00014                                                         :
00015   addRadiatedGluons_( cfg.getParameter<bool>( "addRadiatedGluons" ) ),
00016   src_( cfg.getParameter<edm::InputTag>( "src" ) )
00017 {
00018   // produces a set of gen particle collections following
00019   // the decay branch of top quarks to the first level of 
00020   // stable decay products
00021   produces<reco::GenParticleCollection>(); 
00022   produces<reco::GenParticleCollection>("ME");
00023 }

TopDecaySubset::~TopDecaySubset (  ) 

default destructor

Definition at line 25 of file TopDecaySubset.cc.

00026 {
00027 }


Member Function Documentation

void TopDecaySubset::addDaughters ( int idx,
const reco::GenParticle::const_iterator  part,
reco::GenParticleCollection target,
bool  recursive = true 
) [private]

recursively fill vector for all further dacay particles of a given particle

Definition at line 425 of file TopDecaySubset.cc.

References ptr, and refs_.

Referenced by fromFullListing(), and fromTruncListing().

00426 {
00427   std::vector<int> daughters;
00428   int idxBuffer = idx;
00429   for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
00430     std::auto_ptr<reco::GenParticle> ptr( new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
00431     target.push_back( *ptr );
00432     // increment & push index of daughter
00433     daughters.push_back( ++idx );
00434     // continue recursively if desired
00435     if(recursive){
00436       addDaughters(idx,daughter,target);  
00437     }
00438   }  
00439   if(daughters.size()) {
00440     refs_[ idxBuffer ] = daughters;
00441   }
00442 }

void TopDecaySubset::addRadiation ( int idx,
const reco::GenParticle::const_iterator  part,
reco::GenParticleCollection target 
) [private]

fill vector including all radiations from quarks originating from W/top

Definition at line 405 of file TopDecaySubset.cc.

References ptr, refs_, and TopDecayID::stable.

Referenced by fromFullListing(), and fromTruncListing().

00406 {
00407   std::vector<int> daughters;
00408   int idxBuffer = idx;
00409   for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
00410     if(daughter->status()<=TopDecayID::stable && daughter->pdgId()!=part->pdgId()){
00411       // skip comment lines and make sure that
00412       // the particle is not double counted as 
00413       // dasughter of itself
00414       std::auto_ptr<reco::GenParticle> ptr(  new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
00415       target.push_back( *ptr );
00416       daughters.push_back( ++idx ); //push index of daughter
00417     }
00418   }  
00419   if(daughters.size()) {
00420     refs_[ idxBuffer ] = daughters;
00421   }
00422 }

void TopDecaySubset::clearReferences (  )  [private]

clear references

Definition at line 78 of file TopDecaySubset.cc.

References motherPartIdx_, and refs_.

Referenced by fillOutput().

00079 {
00080   // clear vector of references 
00081   refs_.clear();  
00082   // set idx for mother particles to a start value
00083   // of -1 (the first entry will raise it to 0)
00084   motherPartIdx_=-1;
00085 }

void TopDecaySubset::fillOutput ( edm::Event evt,
const reco::GenParticleCollection src,
FillMode  mode 
) [private]

fill output collection depending on whether the W boson is contained in the generator listing or not

Definition at line 47 of file TopDecaySubset.cc.

References clearReferences(), fillReferences(), fromFullListing(), fromTruncListing(), edm::Event::getRefBeforePut(), kME, kStable, printTarget(), edm::Event::put(), target, TopDecayID::tID, and wInDecayChain().

Referenced by produce().

00048 {
00049   // get ref product from the event
00050   const reco::GenParticleRefProd ref = evt.getRefBeforePut<reco::GenParticleCollection>(); 
00051   // create target vector
00052   std::auto_ptr<reco::GenParticleCollection> target( new reco::GenParticleCollection );
00053 
00054   // clear existing refs
00055   clearReferences();  
00056   // fill output for top branch
00057   wInDecayChain(src, TopDecayID::tID) ? fromFullListing (src, *target, TopDecayID::tID, mode): fromTruncListing(src, *target, TopDecayID::tID, mode);
00058   // fill output for anti-top branch
00059   wInDecayChain(src,-TopDecayID::tID) ? fromFullListing (src, *target,-TopDecayID::tID, mode): fromTruncListing(src, *target,-TopDecayID::tID, mode);
00060   // fill references
00061   fillReferences(ref, *target);
00062   // print full isting of input collection for debuging
00063   // with 'TopDecaySubset_printTarget'
00064   printTarget(*target);
00065 
00066   // write vectors to the event
00067   switch(mode){
00068   case kStable:
00069     evt.put(target);
00070     break;
00071   case kME:
00072     evt.put(target, "ME");
00073     break;
00074   }
00075 }

void TopDecaySubset::fillReferences ( const reco::GenParticleRefProd refProd,
reco::GenParticleCollection target 
) [private]

fill references for output vector

Definition at line 445 of file TopDecaySubset.cc.

References reco::CompositeRefCandidateT< D >::addDaughter(), edm::errors::InvalidReference, p, and refs_.

Referenced by fillOutput().

00446 { 
00447   int idx=0;
00448   for(reco::GenParticleCollection::iterator p=sel.begin(); p!=sel.end(); ++p, ++idx){
00449     //find daughter reference vectors in refs_ and add daughters
00450     std::map<int, std::vector<int> >::const_iterator daughters=refs_.find( idx );
00451     if( daughters!=refs_.end() ){
00452       for(std::vector<int>::const_iterator daughter = daughters->second.begin(); 
00453           daughter!=daughters->second.end(); ++daughter){
00454         reco::GenParticle* part = dynamic_cast<reco::GenParticle* > (&(*p));
00455         if(part == 0){
00456          throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
00457         }
00458         part->addDaughter( reco::GenParticleRef(ref, *daughter) );
00459         sel[*daughter].addMother( reco::GenParticleRef(ref, idx) );
00460       }
00461     }
00462   }
00463 }

void TopDecaySubset::fromFullListing ( const reco::GenParticleCollection src,
reco::GenParticleCollection target,
const int partId,
FillMode  mode 
) [private]

fill output vector with full decay chain with intermediate W bosons

Definition at line 113 of file TopDecaySubset.cc.

References funct::abs(), addDaughters(), addRadiatedGluons_, addRadiation(), TopDecayID::bID, false, TopDecayID::glueID, kME, motherPartIdx_, p4(), refs_, TopDecayID::stable, t, TopDecayID::tauID, getDQMSummary::td, TopDecayID::tID, true, TopDecayID::unfrag, and TopDecayID::WID.

Referenced by fillOutput().

00114 {
00115   int statusFlag;
00116   // determine status flag of the new 
00117   // particle depending on the FillMode
00118   mode == kME ? statusFlag=3 : statusFlag=2;
00119 
00120   for(reco::GenParticleCollection::const_iterator t=src.begin(); t!=src.end(); ++t){
00121     if( t->status() == TopDecayID::unfrag && t->pdgId()==partId ){ 
00122       // if particle is top or anti-top depending 
00123       // on the function's argument partId 
00124       std::auto_ptr<reco::GenParticle> topPtr( new reco::GenParticle( t->threeCharge(), p4( t, statusFlag), t->vertex(), t->pdgId(), statusFlag, false ) );
00125       target.push_back( *topPtr );
00126       ++motherPartIdx_;
00127       // keep the top index for the map to manage the daughter refs
00128       int iTop=motherPartIdx_; 
00129       std::vector<int> topDaughters;
00130       // define the W boson index (to be set later) for the map to 
00131       // manage the daughter refs
00132       int iW = 0;
00133       std::vector<int> wDaughters;
00134       //iterate over top daughters
00135       for(reco::GenParticle::const_iterator td=t->begin(); td!=t->end(); ++td){
00136         if( td->status()==TopDecayID::unfrag && abs( td->pdgId() )<=TopDecayID::bID ){ 
00137           // if particle is beauty or other quark 
00138           std::auto_ptr<reco::GenParticle> bPtr( new reco::GenParticle( td->threeCharge(), p4( td, statusFlag ), td->vertex(), td->pdgId(), statusFlag, false ) );
00139           target.push_back( *bPtr );      
00140           // increment & push index of the top daughter
00141           topDaughters.push_back( ++motherPartIdx_ ); 
00142           if(addRadiatedGluons_){
00143             addRadiation(motherPartIdx_,td,target); 
00144           }
00145         }
00146         if( td->status()==TopDecayID::unfrag && abs( td->pdgId() )==TopDecayID::WID ){ 
00147           // if particle is is W boson
00148           std::auto_ptr<reco::GenParticle> wPtr(  new reco::GenParticle( td->threeCharge(), p4( td, statusFlag), td->vertex(), td->pdgId(), statusFlag, true ) );
00149           target.push_back( *wPtr );
00150           // increment & push index of the top daughter
00151           topDaughters.push_back( ++motherPartIdx_ );
00152           // keep the W idx for the map
00153           iW=motherPartIdx_; 
00154           if(addRadiatedGluons_){
00155             addRadiation(motherPartIdx_,td,target); 
00156           }
00157 
00158           // iterate over W daughters
00159           for(reco::GenParticle::const_iterator wd=td->begin(); wd!=td->end(); ++wd){
00160             // make sure the W daughter is of status unfrag and not the W itself
00161             if( wd->status()==TopDecayID::unfrag && !(abs(wd->pdgId())==TopDecayID::WID) ) {
00162               std::auto_ptr<reco::GenParticle> qPtr( new reco::GenParticle( wd->threeCharge(), p4( wd, statusFlag ), wd->vertex(), wd->pdgId(), statusFlag, false) );
00163               target.push_back( *qPtr );
00164               // increment & push index of the top daughter
00165               wDaughters.push_back( ++motherPartIdx_ );
00166               if(addRadiatedGluons_){
00167                 addRadiation(motherPartIdx_,wd,target); 
00168               }
00169               if( wd->status()==TopDecayID::unfrag && abs( wd->pdgId() )==TopDecayID::tauID ){ 
00170                 // add tau daughters if the particle is a tau pass
00171                 // the daughter of the tau which is of status 2
00172                 addDaughters(motherPartIdx_,wd->begin(),target); 
00173               }
00174             } 
00175           }
00176         }
00177         if(td->status()==TopDecayID::stable && ( td->pdgId()==TopDecayID::glueID || abs(td->pdgId())<TopDecayID::bID)){
00178           // collect additional radiation from the top 
00179           std::auto_ptr<reco::GenParticle> radPtr( new reco::GenParticle( td->threeCharge(), td->p4(), td->vertex(), td->pdgId(), statusFlag, false ) );
00180           target.push_back( *radPtr );  
00181           // increment & push index of the top daughter  
00182           topDaughters.push_back( ++motherPartIdx_ );
00183         }
00184       }
00185       // add potential sisters of the top quark;
00186       // only for top to prevent double counting
00187       if(t->numberOfMothers()>0 && t->pdgId()==TopDecayID::tID){
00188         for(reco::GenParticle::const_iterator ts = t->mother()->begin(); ts!=t->mother()->end(); ++ts){
00189           // loop over all daughters of the top mother i.e.
00190           // the two top quarks and their potential sisters
00191           if(abs(ts->pdgId())!=t->pdgId()){
00192             // add all further particles
00193             // but the two top quarks 
00194             reco::GenParticle* cand = new reco::GenParticle( ts->threeCharge(), ts->p4(), ts->vertex(), ts->pdgId(), ts->status(), false );
00195             std::auto_ptr<reco::GenParticle> sPtr( cand );
00196             target.push_back( *sPtr );
00197             // increment & add 1. generation of daughters
00198             addDaughters(++motherPartIdx_,ts->begin(),target,false);
00199           }
00200         }
00201       }
00202       // fill the map for the administration 
00203       // of daughter indices
00204       refs_[ iTop ] = topDaughters;
00205       refs_[ iW   ] = wDaughters; 
00206     }
00207   }
00208 }

void TopDecaySubset::fromTruncListing ( const reco::GenParticleCollection src,
reco::GenParticleCollection target,
const int partId,
FillMode  mode 
) [private]

fill output vector with full decay chain w/o intermediate W bosons

Definition at line 211 of file TopDecaySubset.cc.

References funct::abs(), addDaughters(), addRadiatedGluons_, addRadiation(), TopDecayID::bID, reco::Particle::charge(), false, TopDecayID::glueID, kME, motherPartIdx_, reco::Particle::p4(), p4(), refs_, TopDecayID::stable, t, TopDecayID::tauID, getDQMSummary::td, reco::Particle::threeCharge(), TopDecayID::tID, true, TopDecayID::unfrag, and reco::Particle::vertex().

Referenced by fillOutput().

00212 {
00213   int statusFlag;
00214   // determine status flag of the new 
00215   // particle depending on the FillMode
00216   mode == kME ? statusFlag=3 : statusFlag=2;
00217 
00218   // needed for W reconstruction from 
00219   // daughters
00220   reco::Particle::Point wVtx;
00221   reco::Particle::Charge wQ=0;
00222   reco::Particle::LorentzVector wP4;
00223   for(reco::GenParticleCollection::const_iterator t=src.begin(); t!=src.end(); ++t){
00224     if( t->status() == TopDecayID::unfrag && t->pdgId()==partId ){ 
00225       // if particle is top or anti-top 
00226       std::auto_ptr<reco::GenParticle> topPtr( new reco::GenParticle( t->threeCharge(), p4( t, statusFlag), t->vertex(), t->pdgId(), statusFlag, false ) );
00227       target.push_back( *topPtr );
00228       ++motherPartIdx_;
00229       // keep the top index for the map to manage the daughter refs
00230       int iTop=motherPartIdx_; 
00231       std::vector<int> topDaughters;
00232       // define the W boson index (to be set later) for the map to 
00233       // manage the daughter refs
00234       int iW = 0;
00235       std::vector<int> wDaughters;
00236       
00237       //iterate over top daughters
00238       for(reco::GenParticle::const_iterator td=t->begin(); td!=t->end(); ++td){
00239         if( td->status()==TopDecayID::unfrag && abs( td->pdgId() )==TopDecayID::bID ){ 
00240           // if particle is beauty or other quark 
00241           std::auto_ptr<reco::GenParticle> bPtr( new reco::GenParticle( td->threeCharge(), p4( td, statusFlag ), td->vertex(), td->pdgId(), statusFlag, false ) );
00242           target.push_back( *bPtr );      
00243           // increment & push index of the top daughter
00244           topDaughters.push_back( ++motherPartIdx_ ); 
00245           if(addRadiatedGluons_){
00246             addRadiation(motherPartIdx_,td,target); 
00247           }
00248         }
00249 
00250         // count all 4-vectors but the b quark 
00251         // as daughters of the W boson
00252         if( td->status()==TopDecayID::unfrag && abs( td->pdgId() )!=TopDecayID::bID ){
00253           reco::GenParticle* cand = new reco::GenParticle( td->threeCharge(), p4( td, statusFlag ), td->vertex(), td->pdgId(), statusFlag, false);
00254           std::auto_ptr<reco::GenParticle> qPtr( cand );
00255           target.push_back( *qPtr );
00256           // increment & push index of the top daughter
00257           wDaughters.push_back( ++motherPartIdx_ );
00258           if(addRadiatedGluons_){
00259             if(abs(td->pdgId())<TopDecayID::tID){
00260               addRadiation(motherPartIdx_,td,target); 
00261             }
00262           }
00263 
00264           // reconstruct the quantities of the W boson from its daughters; take 
00265           // care of the non-trivial charge definition of quarks/leptons and of 
00266           // the non-integer charge of the quarks
00267           if( fabs(td->pdgId() )<TopDecayID::tID ) //quark
00268             wQ += ((td->pdgId()>0)-(td->pdgId()<0))*abs(cand->threeCharge());
00269           else //lepton
00270             wQ +=-((td->pdgId()>0)-(td->pdgId()<0))*abs(cand->charge());
00271           wVtx=cand->vertex();        
00272           wP4+=cand->p4();
00273 
00274           if( td->status()==TopDecayID::unfrag && abs( td->pdgId() )==TopDecayID::tauID ){ 
00275             // add tau daughters if the particle is a tau pass
00276             // the daughter of the tau which is of status 2
00277             addDaughters(motherPartIdx_,td->begin(),target); 
00278           }
00279         }
00280         if( (td+1)==t->end() ){
00281           // the reconstruction of the W boson 
00282           // candidate; now add it to the list
00283           std::auto_ptr<reco::GenParticle> wPtr(  new reco::GenParticle( td->threeCharge(), p4( td, statusFlag), td->vertex(), td->pdgId(), statusFlag, true ) );
00284           target.push_back( *wPtr );
00285           // increment & push index of the top daughter
00286           topDaughters.push_back( ++motherPartIdx_ );
00287           // keep the W idx for the map
00288           iW=motherPartIdx_; 
00289           if(addRadiatedGluons_){
00290             addRadiation(motherPartIdx_,td,target); 
00291           }
00292           
00293           // reset the quantities of the W boson 
00294           // for the next reconstruction cycle
00295           wQ  = 0;
00296           wVtx= reco::Particle::Point(0, 0, 0);
00297           wP4 = reco::Particle::LorentzVector(0, 0, 0, 0);
00298         }
00299         if(td->status()==TopDecayID::stable && ( td->pdgId()==TopDecayID::glueID || abs(td->pdgId())<TopDecayID::bID)){
00300           // collect additional radiation from the top 
00301           std::auto_ptr<reco::GenParticle> radPtr( new reco::GenParticle( td->threeCharge(), td->p4(), td->vertex(), td->pdgId(), statusFlag, false ) );
00302           target.push_back( *radPtr );  
00303           // increment & push index of the top daughter  
00304           topDaughters.push_back( ++motherPartIdx_ );
00305         }
00306       }
00307       // add potential sisters of the top quark;
00308       // only for top to prevent double counting
00309       if(t->numberOfMothers()>0 && t->pdgId()==TopDecayID::tID){
00310         for(reco::GenParticle::const_iterator ts = t->mother()->begin(); ts!=t->mother()->end(); ++ts){
00311           // loop over all daughters of the top mother i.e.
00312           // the two top quarks and their potential sisters
00313           if(abs(ts->pdgId())!=t->pdgId()){
00314             // add all further particles
00315             // but the two top quarks 
00316             reco::GenParticle* cand = new reco::GenParticle( ts->threeCharge(), ts->p4(), ts->vertex(), ts->pdgId(), ts->status(), false );
00317             std::auto_ptr<reco::GenParticle> sPtr( cand );
00318             target.push_back( *sPtr );
00319             // increment & add 1. generation of daughters
00320             addDaughters(++motherPartIdx_,ts->begin(),target,false);
00321           }
00322         }
00323       }
00324       // fill the map for the administration 
00325       // of daughter indices
00326       refs_[ iTop ] = topDaughters;
00327       refs_[ iW   ] = wDaughters; 
00328     }
00329   }
00330 }

reco::Particle::LorentzVector TopDecaySubset::p4 ( const reco::GenParticle::const_iterator  part,
int  statusFlag 
) [private]

calculate lorentz vector from input

Definition at line 374 of file TopDecaySubset.cc.

References funct::abs(), p, p4(), TopDecayID::stable, TopDecayID::unfrag, and TopDecayID::WID.

00375 {
00376   // calculate the four vector for all top daughters from 
00377   // their daughters including additional radiation 
00378   if(statusFlag==TopDecayID::unfrag){
00379     // return 4 momentum as it is
00380     return part->p4();
00381   }
00382   reco::Particle::LorentzVector vec;
00383   for(reco::GenParticle::const_iterator p=part->begin(); p!=part->end(); ++p){
00384     if( p->status()<=TopDecayID::stable && 
00385         abs(p->pdgId ())==TopDecayID::WID){
00386       vec=p->p4();
00387     }
00388     else{
00389       if(p->status()<=TopDecayID::stable){
00390         // sum up the p4 of all stable particles 
00391         // (of status 1 or 2)
00392         vec+=p->p4();
00393       }
00394       else 
00395         if( p->status()==TopDecayID::unfrag)
00396           // if the particle is unfragmented (i.e.
00397           // status 3) decend by one level
00398           vec+=p4(p, statusFlag);   
00399     }
00400   }
00401   return vec;
00402 }

reco::Particle::LorentzVector TopDecaySubset::p4 ( const std::vector< reco::GenParticle >::const_iterator  top,
int  statusFlag 
) [private]

calculate lorentz vector from input (dedicated to top reconstruction)

Definition at line 333 of file TopDecaySubset.cc.

References funct::abs(), p, TopDecayID::unfrag, and TopDecayID::WID.

Referenced by fromFullListing(), fromTruncListing(), and p4().

00334 {
00335   // calculate the four vector for top/anti-top quarks from 
00336   // the W boson and the b quark plain or including all 
00337   // additional radiation depending on switch 'plain'
00338   if(statusFlag==TopDecayID::unfrag){
00339     // return 4 momentum as it is
00340     return top->p4();
00341   }
00342   reco::Particle::LorentzVector vec;
00343   for(reco::GenParticle::const_iterator p=top->begin(); p!=top->end(); ++p){
00344     if( p->status() == TopDecayID::unfrag ){
00345       // decend by one level for each
00346       // status 3 particle on the way
00347       vec+=p4( p, statusFlag );
00348     }
00349     else{ 
00350       if( abs(top->pdgId())==TopDecayID::WID ){
00351         // in case of a W boson skip the status 
00352         // 2 particle to prevent double counting
00353         if( abs(p->pdgId())!=TopDecayID::WID ) 
00354           vec+=p->p4();
00355       } 
00356       else{
00357         // add all four vectors for each stable
00358         // particle (status 1 or 2) on the way 
00359         // else
00360         vec+=p->p4();
00361         if( vec.mass()-top->mass()>0 ){
00362           // continue adding up gluons and qqbar pairs on the top 
00363           // line until the nominal top mass is reached; then break 
00364           // in order to prevent picking up virtualities
00365           break;
00366         }
00367       }
00368     }
00369   }
00370   return vec;
00371 }

void TopDecaySubset::printSource ( const reco::GenParticleCollection src  )  [private]

print the whole listing if particle with pdgId is contained in the top decay chain

Definition at line 466 of file TopDecaySubset.cc.

References d, kMAX, funct::log(), p, t, and TopDecayID::tID.

Referenced by produce().

00467 {
00468   edm::LogVerbatim log("TopDecaySubset_printSource");
00469   log << "\n   idx   pdg   stat      px          py         pz             mass          daughter pdg's  "
00470       << "\n===========================================================================================\n";
00471 
00472   for(unsigned int t=0; t<src.size(); ++t){
00473     if( src[t].pdgId()==TopDecayID::tID ){
00474       // restrict to the top in order 
00475       // to have it shown only once 
00476       int idx=0;
00477       for(reco::GenParticleCollection::const_iterator p=src.begin(); p!=src.end(); ++p, ++idx){
00478         // loop the top daughters
00479         log << std::right << std::setw( 5) << idx
00480             << std::right << std::setw( 7) << src[idx].pdgId()
00481             << std::right << std::setw( 5) << src[idx].status() << "  "
00482             << std::right << std::setw(10) << std::setprecision( 6 ) << src[idx].p4().x() << "  "       
00483             << std::right << std::setw(10) << std::setprecision( 6 ) << src[idx].p4().y() << "  "       
00484             << std::right << std::setw(10) << std::setprecision( 6 ) << src[idx].p4().z() << "  "       
00485             << std::right << std::setw(15) << std::setprecision( 6 ) << src[idx].p4().mass() 
00486             << "   ";
00487         // search for potential daughters; if they exits 
00488         // print the daughter to the screen in the last 
00489         // column of the table separated by ','
00490         TString pdgIds;
00491         unsigned int jdx=0;
00492         for(reco::GenParticle::const_iterator d=p->begin(); d!=p->end(); ++d, ++jdx){
00493           if(jdx<kMAX){
00494             pdgIds+=d->pdgId();
00495             if(d+1 != p->end()){
00496               pdgIds+= ",";
00497             }
00498           }
00499           else{
00500             pdgIds+="...(";
00501             pdgIds+= p->numberOfDaughters();
00502             pdgIds+=")";
00503             break;
00504           }
00505         }
00506         if(idx>0){
00507           log << std::setfill( ' ' ) << std::right << std::setw(15) << pdgIds; 
00508           log << "\n";
00509         }
00510         else{
00511           log << std::setfill( ' ' ) << std::right << std::setw(15) << "-\n";
00512         }
00513       }
00514     }
00515   }
00516 }

void TopDecaySubset::printTarget ( reco::GenParticleCollection target  )  [private]

print the whole decay chain if particle with pdgId is contained in the top decay chain

Definition at line 519 of file TopDecaySubset.cc.

References d, kMAX, funct::log(), p, refs_, t, and TopDecayID::tID.

Referenced by fillOutput().

00520 {
00521   edm::LogVerbatim log("TopDecaySubset_printTarget");
00522   log << "\n   idx   pdg   stat      px          py         pz             mass          daughter pdg's  "
00523       << "\n===========================================================================================\n";
00524 
00525   for(unsigned int t=0; t<sel.size(); ++t){
00526     if( sel[t].pdgId()==TopDecayID::tID ){
00527       // restrict to the top in order 
00528       // to have it shown only once      
00529       int idx=0;
00530       for(reco::GenParticleCollection::iterator p=sel.begin(); p!=sel.end(); ++p, ++idx){
00531         // loop the top daughters
00532         log << std::right << std::setw( 5) << idx
00533             << std::right << std::setw( 7) << sel[idx].pdgId()
00534             << std::right << std::setw( 5) << sel[idx].status() << "  "
00535             << std::right << std::setw(10) << std::setprecision( 6 ) << sel[idx].p4().x() << "  "       
00536             << std::right << std::setw(10) << std::setprecision( 6 ) << sel[idx].p4().y() << "  "       
00537             << std::right << std::setw(10) << std::setprecision( 6 ) << sel[idx].p4().z() << "  "       
00538             << std::right << std::setw(15) << std::setprecision( 6 ) << sel[idx].p4().mass() 
00539             << "   ";
00540         // search for potential daughters; if they exits 
00541         // print the daughter to the screen in the last 
00542         // column of the table separated by ','
00543         TString pdgIds;
00544         std::map<int, std::vector<int> >::const_iterator daughters=refs_.find( idx );
00545         if( daughters!=refs_.end() ){
00546           unsigned int jdx=0;
00547           for(std::vector<int>::const_iterator d = daughters->second.begin(); d!=daughters->second.end(); ++d, ++jdx){
00548             if(jdx<kMAX){
00549               pdgIds+=sel[*d].pdgId();
00550               if(d+1 != daughters->second.end()){
00551                 pdgIds+= ",";
00552               }
00553             }
00554             else{
00555               pdgIds+="...(";
00556               pdgIds+= daughters->second.size();
00557               pdgIds+=")";
00558               break;
00559             }
00560           }
00561           log << std::setfill( ' ' ) << std::right << std::setw(15) << pdgIds; 
00562           log << "\n";
00563         }
00564         else{
00565           log << std::setfill( ' ' ) << std::right << std::setw(15) << "-\n";
00566         }
00567       }
00568     }
00569   }
00570 }

void TopDecaySubset::produce ( edm::Event evt,
const edm::EventSetup setup 
) [virtual]

write TopDecaySubset into the event

Implements edm::EDProducer.

Definition at line 30 of file TopDecaySubset.cc.

References fillOutput(), edm::Event::getByLabel(), kME, kStable, printSource(), HLT_VtxMuL3::src, and src_.

00031 {     
00032   // get source collection
00033   edm::Handle<reco::GenParticleCollection> src;
00034   evt.getByLabel(src_, src);
00035 
00036   // print full listing of input collection for debuging
00037   // with 'TopDecaySubset_printSource'
00038   printSource(*src);
00039 
00040   // fill output vectors with references
00041   fillOutput (evt, *src, kME      );
00042   fillOutput (evt, *src, kStable  );
00043 }

bool TopDecaySubset::wInDecayChain ( const reco::GenParticleCollection src,
const int partId 
) [private]

check whether the W boson is contained in the original gen particle listing or not

Definition at line 88 of file TopDecaySubset.cc.

References funct::abs(), lhef::isContained(), t, getDQMSummary::td, TopDecayID::unfrag, and TopDecayID::WID.

Referenced by fillOutput().

00089 {
00090   bool isContained=false;
00091   for(reco::GenParticleCollection::const_iterator t=src.begin(); t!=src.end(); ++t){
00092     if( t->status() == TopDecayID::unfrag && t->pdgId()==partId ){ 
00093       for(reco::GenParticle::const_iterator td=t->begin(); td!=t->end(); ++td){
00094         if( td->status()==TopDecayID::unfrag && abs( td->pdgId() )==TopDecayID::WID ){ 
00095           isContained=true;
00096           break;
00097         }
00098       }
00099     }
00100   }
00101   if( !isContained ){
00102     edm::LogWarning( "decayChain" )
00103       << " W boson is not contained in decay chain in the gen particle listing \n"
00104       << " of the input collection. A status 2 equivalent W candidate will be  \n"
00105       << " re-reconstructed from the W daughters but the hadronization of the  \n"
00106       << " W boson might be screwed up. Contact an expert for the generator in \n"
00107       << " use to assure that what you are doing makes sense.";     
00108   }
00109   return isContained;
00110 }


Member Data Documentation

bool TopDecaySubset::addRadiatedGluons_ [private]

add radiated gluons or not?

Definition at line 76 of file TopDecaySubset.h.

Referenced by fromFullListing(), and fromTruncListing().

int TopDecaySubset::motherPartIdx_ [private]

index in new evt listing of parts with daughters; has to be set to -1 in produce to deliver consistent results!

Definition at line 72 of file TopDecaySubset.h.

Referenced by clearReferences(), fromFullListing(), and fromTruncListing().

std::map<int,std::vector<int> > TopDecaySubset::refs_ [private]

management of daughter indices for fillRefs

Definition at line 74 of file TopDecaySubset.h.

Referenced by addDaughters(), addRadiation(), clearReferences(), fillReferences(), fromFullListing(), fromTruncListing(), and printTarget().

edm::InputTag TopDecaySubset::src_ [private]

input tag for the genParticle source

Definition at line 78 of file TopDecaySubset.h.

Referenced by produce().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:33:41 2009 for CMSSW by  doxygen 1.5.4