00001 #include "TString.h"
00002 #include "FWCore/Utilities/interface/EDMException.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004
00005 #include "AnalysisDataFormats/TopObjects/interface/TopGenEvent.h"
00006 #include "TopQuarkAnalysis/TopEventProducers/interface/TopDecaySubset.h"
00007
00009 TopDecaySubset::TopDecaySubset(const edm::ParameterSet& cfg):
00010 src_ (cfg.getParameter<edm::InputTag>("src")),
00011 addRadiation_(cfg.getParameter<bool>("addRadiation")),
00012 showerModel_ (kStart)
00013 {
00014
00015
00016 if( cfg.getParameter<std::string>( "fillMode" ) == "kME" ){ fillMode_=kME; }
00017 if( cfg.getParameter<std::string>( "fillMode" ) == "kStable" ){ fillMode_=kStable; }
00018
00019
00020
00021 produces<reco::GenParticleCollection>();
00022 }
00023
00025 TopDecaySubset::~TopDecaySubset()
00026 {
00027 }
00028
00030 void
00031 TopDecaySubset::produce(edm::Event& event, const edm::EventSetup& setup)
00032 {
00033
00034 edm::Handle<reco::GenParticleCollection> src;
00035 event.getByLabel(src_, src);
00036
00037
00038 std::vector<const reco::GenParticle*> tops = findTops(*src);
00039
00040
00041 if(showerModel_==kStart)
00042 showerModel_=checkShowerModel(tops);
00043
00044
00045 checkWBosons(tops);
00046
00047
00048 std::auto_ptr<reco::GenParticleCollection> target( new reco::GenParticleCollection );
00049
00050 const reco::GenParticleRefProd ref = event.getRefBeforePut<reco::GenParticleCollection>();
00051
00052 clearReferences();
00053
00054 fillListing(tops, *target);
00055
00056 fillReferences(ref, *target);
00057
00058
00059 event.put(target);
00060 }
00061
00063 std::vector<const reco::GenParticle*>
00064 TopDecaySubset::findTops(const reco::GenParticleCollection& parts)
00065 {
00066 std::vector<const reco::GenParticle*> tops;
00067 for(reco::GenParticleCollection::const_iterator t=parts.begin(); t!=parts.end(); ++t){
00068 if( std::abs(t->pdgId())==TopDecayID::tID && t->status()==TopDecayID::unfrag )
00069 tops.push_back( &(*t) );
00070 }
00071 return tops;
00072 }
00073
00075 TopDecaySubset::ShowerModel
00076 TopDecaySubset::checkShowerModel(const std::vector<const reco::GenParticle*>& tops) const
00077 {
00078 for(std::vector<const reco::GenParticle*>::const_iterator it=tops.begin(); it!=tops.end(); ++it){
00079 const reco::GenParticle* top = *it;
00080
00081
00082
00083 if( top->numberOfDaughters()==1){
00084 if( top->begin()->pdgId()==top->pdgId() && top->begin()->status()==TopDecayID::stable && top->begin()->numberOfDaughters()>1)
00085 return kHerwig;
00086 }
00087
00088
00089
00090 if( top->numberOfDaughters()>1 ){
00091 bool containsWBoson=false, containsQuarkDaughter=false;
00092 for(reco::GenParticle::const_iterator td=top->begin(); td!=top->end(); ++td){
00093 if( std::abs(td->pdgId ()) <TopDecayID::tID ) containsQuarkDaughter=true;
00094 if( std::abs(td->pdgId ())==TopDecayID::WID ) containsWBoson=true;
00095 }
00096 if(containsQuarkDaughter && containsWBoson)
00097 return kPythia;
00098 }
00099 }
00100
00101 if(tops.size()==0)
00102 edm::LogWarning("decayChain")
00103 << " Failed to find top quarks in decay chain. Will assume that this a \n"
00104 << " non-top sample and produce an empty decaySubset.\n";
00105 else
00106 throw edm::Exception(edm::errors::LogicError,
00107 " Can not find back any of the supported hadronization models. Models \n"
00108 " which are supported are: \n"
00109 " Pythia LO(+PS): Top(status 3) --> WBoson(status 3), Quark(status 3)\n"
00110 " Herwig NLO(+PS): Top(status 2) --> Top(status 3) --> Top(status 2) \n");
00111 return kNone;
00112 }
00113
00115 void
00116 TopDecaySubset::checkWBosons(std::vector<const reco::GenParticle*>& tops) const
00117 {
00118 unsigned nTops = tops.size();
00119 for(std::vector<const reco::GenParticle*>::iterator it=tops.begin(); it!=tops.end();) {
00120 const reco::GenParticle* top = *it;
00121 bool isContained=false;
00122 bool expectedStatus=false;
00123 for(reco::GenParticle::const_iterator td=((showerModel_==kPythia) ? top->begin() : top->begin()->begin());
00124 td!=((showerModel_==kPythia) ? top->end() : top->begin()->end());
00125 ++td){
00126 if(std::abs(td->pdgId())==TopDecayID::WID) {
00127 isContained=true;
00128 if( ((showerModel_==kPythia) ? td->status()==TopDecayID::unfrag : td->status()==TopDecayID::stable) ) {
00129 expectedStatus=true;
00130 break;
00131 }
00132 }
00133 }
00134 if(!expectedStatus) {
00135 it=tops.erase(it);
00136 if(isContained)
00137 edm::LogWarning("decayChain")
00138 << " W boson does not have the expected status. This happens, e.g., \n"
00139 << " with MC@NLO in the case of additional ttbar pairs from radiated \n"
00140 << " gluons. We hope everything is fine, remove the correspdonding top \n"
00141 << " quark from our list since it is not part of the primary ttbar pair \n"
00142 << " and try to continue. \n";
00143 }
00144 else
00145 it++;
00146 }
00147 if(tops.size()==0 && nTops!=0)
00148 throw edm::Exception(edm::errors::LogicError,
00149 " Did not find a W boson with appropriate status for any of the top \n"
00150 " quarks in this event. This means that the hadronization of the W \n"
00151 " boson might be screwed up or there is another problem with the \n"
00152 " particle listing in this MC sample. Please contact an expert. \n");
00153 }
00154
00156 void
00157 TopDecaySubset::fillListing(const std::vector<const reco::GenParticle*>& tops, reco::GenParticleCollection& target)
00158 {
00159 unsigned int statusFlag;
00160
00161
00162 fillMode_ == kME ? statusFlag=3 : statusFlag=2;
00163
00164 for(std::vector<const reco::GenParticle*>::const_iterator it=tops.begin(); it!=tops.end(); ++it){
00165 const reco::GenParticle* t = *it;
00166
00167 std::auto_ptr<reco::GenParticle> topPtr( new reco::GenParticle( t->threeCharge(), p4(it, statusFlag), t->vertex(), t->pdgId(), statusFlag, false ) );
00168 target.push_back( *topPtr );
00169 ++motherPartIdx_;
00170
00171 int iTop=motherPartIdx_;
00172 std::vector<int> topDaughters;
00173
00174
00175 int iW = 0;
00176 std::vector<int> wDaughters;
00177
00178 for(reco::GenParticle::const_iterator td=((showerModel_==kPythia)?t->begin():t->begin()->begin()); td!=((showerModel_==kPythia)?t->end():t->begin()->end()); ++td){
00179 if( td->status()==TopDecayID::unfrag && std::abs( td->pdgId() )<=TopDecayID::bID ){
00180
00181 std::auto_ptr<reco::GenParticle> bPtr( new reco::GenParticle( td->threeCharge(), p4( td, statusFlag ), td->vertex(), td->pdgId(), statusFlag, false ) );
00182 target.push_back( *bPtr );
00183
00184 topDaughters.push_back( ++motherPartIdx_ );
00185 if(addRadiation_){
00186 addRadiation(motherPartIdx_,td,target);
00187 }
00188 }
00189 reco::GenParticle::const_iterator buffer = (showerModel_==kPythia)?td:td->begin();
00190 if( buffer->status()==TopDecayID::unfrag && std::abs( buffer->pdgId() )==TopDecayID::WID ){
00191
00192 std::auto_ptr<reco::GenParticle> wPtr( new reco::GenParticle( buffer->threeCharge(), p4( buffer, statusFlag), buffer->vertex(), buffer->pdgId(), statusFlag, true ) );
00193 target.push_back( *wPtr );
00194
00195 topDaughters.push_back( ++motherPartIdx_ );
00196
00197 iW=motherPartIdx_;
00198 if(addRadiation_){
00199 addRadiation(motherPartIdx_,buffer,target);
00200 }
00201
00202 for(reco::GenParticle::const_iterator wd=((showerModel_==kPythia)?buffer->begin():buffer->begin()->begin()); wd!=((showerModel_==kPythia)?buffer->end():buffer->begin()->end()); ++wd){
00203
00204 if( wd->status()==TopDecayID::unfrag && !(std::abs(wd->pdgId())==TopDecayID::WID) ) {
00205 std::auto_ptr<reco::GenParticle> qPtr( new reco::GenParticle( wd->threeCharge(), p4( wd, statusFlag ), wd->vertex(), wd->pdgId(), statusFlag, false) );
00206 target.push_back( *qPtr );
00207
00208 wDaughters.push_back( ++motherPartIdx_ );
00209 if( wd->status()==TopDecayID::unfrag && std::abs( wd->pdgId() )==TopDecayID::tauID ){
00210
00211
00212
00213
00214
00215
00216 addDaughters(motherPartIdx_,wd,target);
00217 }
00218 }
00219 }
00220 }
00221 if(addRadiation_ && buffer->status()==TopDecayID::stable && ( buffer->pdgId()==TopDecayID::glueID || std::abs(buffer->pdgId())<TopDecayID::bID)){
00222
00223 std::auto_ptr<reco::GenParticle> radPtr( new reco::GenParticle( buffer->threeCharge(), buffer->p4(), buffer->vertex(), buffer->pdgId(), statusFlag, false ) );
00224 target.push_back( *radPtr );
00225
00226 topDaughters.push_back( ++motherPartIdx_ );
00227 }
00228 }
00229
00230
00231 if(t->numberOfMothers()>0 && t->pdgId()==TopDecayID::tID){
00232 for(reco::GenParticle::const_iterator ts = t->mother()->begin(); ts!=t->mother()->end(); ++ts){
00233
00234
00235 if(std::abs(ts->pdgId())!=t->pdgId() && ts->pdgId()!=t->mother()->pdgId()){
00236
00237
00238 reco::GenParticle* cand = new reco::GenParticle( ts->threeCharge(), ts->p4(), ts->vertex(), ts->pdgId(), ts->status(), false );
00239 std::auto_ptr<reco::GenParticle> sPtr( cand );
00240 target.push_back( *sPtr );
00241 if(ts->begin()!=ts->end()){
00242
00243
00244 addDaughters(++motherPartIdx_,ts->begin(),target,false);
00245 }
00246 }
00247 }
00248 }
00249
00250
00251 refs_[ iTop ] = topDaughters;
00252 refs_[ iW ] = wDaughters;
00253 }
00254 }
00255
00257 reco::Particle::LorentzVector
00258 TopDecaySubset::p4(const std::vector<const reco::GenParticle*>::const_iterator top, int statusFlag)
00259 {
00260
00261
00262
00263 if(statusFlag==TopDecayID::unfrag){
00264
00265 return (*top)->p4();
00266 }
00267 reco::Particle::LorentzVector vec;
00268 for(reco::GenParticle::const_iterator p=(*top)->begin(); p!=(*top)->end(); ++p){
00269 if( p->status() == TopDecayID::unfrag ){
00270
00271
00272 vec+=p4( p, statusFlag );
00273 }
00274 else{
00275 if( std::abs((*top)->pdgId())==TopDecayID::WID ){
00276
00277
00278 if( std::abs(p->pdgId())!=TopDecayID::WID )
00279 vec+=p->p4();
00280 }
00281 else{
00282
00283
00284 vec+=p->p4();
00285 if( vec.mass()-(*top)->mass()>0 ){
00286
00287
00288
00289 break;
00290 }
00291 }
00292 }
00293 }
00294 return vec;
00295 }
00296
00298 reco::Particle::LorentzVector
00299 TopDecaySubset::p4(const reco::GenParticle::const_iterator part, int statusFlag)
00300 {
00301
00302
00303 if(statusFlag==TopDecayID::unfrag){
00304
00305 return part->p4();
00306 }
00307 reco::Particle::LorentzVector vec;
00308 for(reco::GenParticle::const_iterator p=part->begin(); p!=part->end(); ++p){
00309 if( p->status()<=TopDecayID::stable && std::abs(p->pdgId ())==TopDecayID::WID){
00310 vec=p->p4();
00311 }
00312 else{
00313 if(p->status()<=TopDecayID::stable){
00314
00315
00316 vec+=p->p4();
00317 }
00318 else{
00319 if( p->status()==TopDecayID::unfrag){
00320
00321
00322 vec+=p4(p, statusFlag);
00323 }
00324 }
00325 }
00326 }
00327 return vec;
00328 }
00329
00331 void
00332 TopDecaySubset::addRadiation(int& idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection& target)
00333 {
00334 std::vector<int> daughters;
00335 int idxBuffer = idx;
00336 for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
00337 if(daughter->status()<=TopDecayID::stable && daughter->pdgId()!=part->pdgId()){
00338
00339
00340
00341 std::auto_ptr<reco::GenParticle> ptr( new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
00342 target.push_back( *ptr );
00343 daughters.push_back( ++idx );
00344 }
00345 }
00346 if(daughters.size()) {
00347 refs_[ idxBuffer ] = daughters;
00348 }
00349 }
00350
00352 void
00353 TopDecaySubset::addDaughters(int& idx, const reco::GenParticle::const_iterator part, reco::GenParticleCollection& target, bool recursive)
00354 {
00355 std::vector<int> daughters;
00356 int idxBuffer = idx;
00357 for(reco::GenParticle::const_iterator daughter=part->begin(); daughter!=part->end(); ++daughter){
00358 std::auto_ptr<reco::GenParticle> ptr( new reco::GenParticle( daughter->threeCharge(), daughter->p4(), daughter->vertex(), daughter->pdgId(), daughter->status(), false) );
00359 target.push_back( *ptr );
00360
00361 daughters.push_back( ++idx );
00362
00363 if(recursive){
00364 addDaughters(idx,daughter,target);
00365 }
00366 }
00367 if(daughters.size()) {
00368 refs_[ idxBuffer ] = daughters;
00369 }
00370 }
00371
00373 void
00374 TopDecaySubset::clearReferences()
00375 {
00376
00377 refs_.clear();
00378
00379
00380 motherPartIdx_=-1;
00381 }
00382
00384 void
00385 TopDecaySubset::fillReferences(const reco::GenParticleRefProd& ref, reco::GenParticleCollection& sel)
00386 {
00387 int idx=0;
00388 for(reco::GenParticleCollection::iterator p=sel.begin(); p!=sel.end(); ++p, ++idx){
00389
00390 std::map<int, std::vector<int> >::const_iterator daughters=refs_.find( idx );
00391 if( daughters!=refs_.end() ){
00392 for(std::vector<int>::const_iterator daughter = daughters->second.begin();
00393 daughter!=daughters->second.end(); ++daughter){
00394 reco::GenParticle* part = dynamic_cast<reco::GenParticle* > (&(*p));
00395 if(part == 0){
00396 throw edm::Exception( edm::errors::InvalidReference, "Not a GenParticle" );
00397 }
00398 part->addDaughter( reco::GenParticleRef(ref, *daughter) );
00399 sel[*daughter].addMother( reco::GenParticleRef(ref, idx) );
00400 }
00401 }
00402 }
00403 }