CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/TopQuarkAnalysis/TopJetCombination/src/TtSemiLepHypothesis.cc

Go to the documentation of this file.
00001 #include "CommonTools/CandUtils/interface/AddFourMomenta.h"
00002 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00003 #include "DataFormats/PatCandidates/interface/Electron.h"
00004 #include "DataFormats/PatCandidates/interface/Muon.h"
00005 #include "TopQuarkAnalysis/TopTools/interface/MEzCalculator.h"
00006 
00007 #include "TopQuarkAnalysis/TopJetCombination/interface/TtSemiLepHypothesis.h"
00008 
00010 TtSemiLepHypothesis::TtSemiLepHypothesis(const edm::ParameterSet& cfg):
00011   jets_(cfg.getParameter<edm::InputTag>("jets")),
00012   leps_(cfg.getParameter<edm::InputTag>("leps")),
00013   mets_(cfg.getParameter<edm::InputTag>("mets")),
00014   nJetsConsidered_(cfg.getParameter<edm::InputTag>("nJetsConsidered")),
00015   numberOfRealNeutrinoSolutions_(-1),
00016   lightQ_(0), lightQBar_(0), hadronicB_(0), 
00017   leptonicB_(0), neutrino_(0), lepton_(0)
00018 {
00019   getMatch_ = false;
00020   if( cfg.exists("match") ) {
00021     getMatch_ = true;
00022     match_ = cfg.getParameter<edm::InputTag>("match");
00023   }
00024   if( cfg.exists("neutrinoSolutionType") )
00025     neutrinoSolutionType_ = cfg.getParameter<int>("neutrinoSolutionType");
00026   else
00027     neutrinoSolutionType_ = -1;
00028   if( cfg.exists("jetCorrectionLevel") ) {
00029     jetCorrectionLevel_ = cfg.getParameter<std::string>("jetCorrectionLevel");
00030   }
00031   produces<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >();
00032   produces<int>("Key");
00033   produces<int>("NumberOfRealNeutrinoSolutions");
00034   produces<int>("NumberOfConsideredJets");
00035 }
00036 
00038 TtSemiLepHypothesis::~TtSemiLepHypothesis()
00039 {
00040   if( lightQ_    ) delete lightQ_;
00041   if( lightQBar_ ) delete lightQBar_;
00042   if( hadronicB_ ) delete hadronicB_;
00043   if( leptonicB_ ) delete leptonicB_;
00044   if( neutrino_  ) delete neutrino_;
00045   if( lepton_    ) delete lepton_;
00046 }
00047 
00049 void
00050 TtSemiLepHypothesis::produce(edm::Event& evt, const edm::EventSetup& setup)
00051 {
00052   edm::Handle<std::vector<pat::Jet> > jets;
00053   evt.getByLabel(jets_, jets);
00054   
00055   edm::Handle<edm::View<reco::RecoCandidate> > leps;
00056   evt.getByLabel(leps_, leps);
00057 
00058   edm::Handle<std::vector<pat::MET> > mets;
00059   evt.getByLabel(mets_, mets);
00060 
00061   edm::Handle<int> nJetsConsidered;
00062   evt.getByLabel(nJetsConsidered_, nJetsConsidered);
00063 
00064   std::vector<std::vector<int> > matchVec;
00065   if( getMatch_ ) {
00066     edm::Handle<std::vector<std::vector<int> > > matchHandle;
00067     evt.getByLabel(match_, matchHandle);
00068     matchVec = *matchHandle;
00069   }
00070   else {
00071     std::vector<int> dummyMatch;
00072     for(unsigned int i = 0; i < 4; ++i) 
00073       dummyMatch.push_back( -1 );
00074     matchVec.push_back( dummyMatch );
00075   }
00076 
00077   // declare auto_ptr for products
00078   std::auto_ptr<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >
00079     pOut( new std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > );
00080   std::auto_ptr<int> pKey(new int);
00081   std::auto_ptr<int> pNeutrinoSolutions(new int);
00082   std::auto_ptr<int> pJetsConsidered(new int);
00083 
00084   // go through given vector of jet combinations
00085   unsigned int idMatch = 0;
00086   typedef std::vector<std::vector<int> >::iterator MatchVecIterator;
00087   for(MatchVecIterator match = matchVec.begin(); match != matchVec.end(); ++match) {
00088     // reset pointers
00089     resetCandidates();
00090     // build hypothesis
00091     buildHypo(evt, leps, mets, jets, *match, idMatch++);
00092     pOut->push_back( std::make_pair(hypo(), *match) );
00093   }
00094   // feed out hyps and matches
00095   evt.put(pOut);
00096 
00097   // build and feed out key
00098   buildKey();
00099   *pKey=key();
00100   evt.put(pKey, "Key");
00101 
00102   // feed out number of real neutrino solutions
00103   *pNeutrinoSolutions=numberOfRealNeutrinoSolutions_;
00104   evt.put(pNeutrinoSolutions, "NumberOfRealNeutrinoSolutions");
00105 
00106   // feed out number of considered jets
00107   *pJetsConsidered=*nJetsConsidered;
00108   evt.put(pJetsConsidered, "NumberOfConsideredJets");
00109 }
00110 
00112 void
00113 TtSemiLepHypothesis::resetCandidates()
00114 {
00115   numberOfRealNeutrinoSolutions_ = -1;
00116   lightQ_    = 0;
00117   lightQBar_ = 0;
00118   hadronicB_ = 0;
00119   leptonicB_ = 0;
00120   neutrino_  = 0;
00121   lepton_    = 0;
00122 }
00123 
00125 reco::CompositeCandidate
00126 TtSemiLepHypothesis::hypo()
00127 {
00128   // check for sanity of the hypothesis
00129   if( !lightQ_ || !lightQBar_ || !hadronicB_ || 
00130       !leptonicB_ || !neutrino_ || !lepton_ )
00131     return reco::CompositeCandidate();
00132   
00133   // setup transient references
00134   reco::CompositeCandidate hyp, hadTop, hadW, lepTop, lepW;
00135 
00136   AddFourMomenta addFourMomenta;  
00137   // build up the top branch that decays leptonically
00138   lepW  .addDaughter(*lepton_,   TtSemiLepDaughter::Lep    );
00139   lepW  .addDaughter(*neutrino_, TtSemiLepDaughter::Nu     );
00140   addFourMomenta.set( lepW );
00141   lepTop.addDaughter( lepW,      TtSemiLepDaughter::LepW   );
00142   lepTop.addDaughter(*leptonicB_,TtSemiLepDaughter::LepB   );
00143   addFourMomenta.set( lepTop );
00144   
00145   // build up the top branch that decays hadronically
00146   hadW  .addDaughter(*lightQ_,   TtSemiLepDaughter::HadP   );
00147   hadW  .addDaughter(*lightQBar_,TtSemiLepDaughter::HadQ   );
00148   addFourMomenta.set( hadW );
00149   hadTop.addDaughter( hadW,      TtSemiLepDaughter::HadW   );
00150   hadTop.addDaughter(*hadronicB_,TtSemiLepDaughter::HadB   );
00151   addFourMomenta.set( hadTop );
00152 
00153   // build ttbar hypotheses
00154   hyp.addDaughter( lepTop,       TtSemiLepDaughter::LepTop );
00155   hyp.addDaughter( hadTop,       TtSemiLepDaughter::HadTop );
00156   addFourMomenta.set( hyp );
00157 
00158   return hyp;
00159 }
00160 
00162 WDecay::LeptonType
00163 TtSemiLepHypothesis::leptonType(const reco::RecoCandidate* cand)
00164 {
00165   // check whetherwe are dealing with a reco muon or a reco electron
00166   WDecay::LeptonType type = WDecay::kNone;
00167   if( dynamic_cast<const reco::Muon*>(cand) ){
00168     type = WDecay::kMuon;
00169   } 
00170   else if( dynamic_cast<const reco::GsfElectron*>(cand) ){
00171     type = WDecay::kElec;
00172   }
00173   return type;
00174 }
00175 
00177 std::string
00178 TtSemiLepHypothesis::jetCorrectionLevel(const std::string& quarkType)
00179 {
00180   // jetCorrectionLevel was not configured
00181   if(jetCorrectionLevel_.empty())
00182     throw cms::Exception("Configuration")
00183       << "Unconfigured jetCorrectionLevel. Please use an appropriate, non-empty string.\n";
00184 
00185   // quarkType is unknown
00186   if( !(quarkType=="wQuarkMix" ||
00187         quarkType=="udsQuark" ||
00188         quarkType=="cQuark" ||
00189         quarkType=="bQuark") )
00190     throw cms::Exception("Configuration")
00191       << quarkType << " is unknown as a quarkType for the jetCorrectionLevel.\n";
00192 
00193   // combine correction level; start with a ':' even if 
00194   // there is no flavor tag to be added, as it is needed
00195   // by setCandidate to disentangle the correction tag 
00196   // from a potential flavor tag, which can be empty
00197   std::string level=jetCorrectionLevel_+":";
00198   if( level=="L5Flavor:" || level=="L6UE:" || level=="L7Parton:" ){
00199     if(quarkType=="wQuarkMix"){level+="wMix";}
00200     if(quarkType=="udsQuark" ){level+="uds";}
00201     if(quarkType=="cQuark"   ){level+="charm";}
00202     if(quarkType=="bQuark"   ){level+="bottom";}
00203   }
00204   else{
00205     level+="none";
00206   }
00207   return level;
00208 }
00209 
00211 void 
00212 TtSemiLepHypothesis::setCandidate(const edm::Handle<std::vector<pat::Jet> >& handle, const int& idx,
00213                                   reco::ShallowClonePtrCandidate*& clone, const std::string& correctionLevel)
00214 {
00215   edm::Ptr<pat::Jet> ptr = edm::Ptr<pat::Jet>(handle, idx);
00216   // disentangle the correction from the potential flavor tag 
00217   // by the separating ':'; the flavor tag can be empty though
00218   std::string step   = correctionLevel.substr(0,correctionLevel.find(":"));
00219   std::string flavor = correctionLevel.substr(1+correctionLevel.find(":"));
00220   float corrFactor = 1.;
00221   if(flavor=="wMix")
00222     corrFactor = 0.75*ptr->jecFactor(step, "uds") + 0.25*ptr->jecFactor(step, "charm");
00223   else
00224     corrFactor = ptr->jecFactor(step, flavor);
00225   clone = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), ptr->p4()*corrFactor, ptr->vertex() );
00226 }
00227 
00229 void TtSemiLepHypothesis::setNeutrino(const edm::Handle<std::vector<pat::MET> >& met,
00230                                       const edm::Handle<edm::View<reco::RecoCandidate> >& leps, const int& idx, const int& type)
00231 {
00232   edm::Ptr<pat::MET> ptr = edm::Ptr<pat::MET>(met, 0);
00233   MEzCalculator mez;
00234   mez.SetMET( *(met->begin()) );
00235   if( leptonType( &(leps->front()) ) == WDecay::kMuon )
00236     mez.SetLepton( (*leps)[idx], true );
00237   else if( leptonType( &(leps->front()) ) == WDecay::kElec )
00238     mez.SetLepton( (*leps)[idx], false );
00239   else
00240     throw cms::Exception("UnimplementedFeature") << "Type of lepton given together with MET for solving neutrino kinematics is neither muon nor electron.\n";
00241   double pz = mez.Calculate(type);
00242   numberOfRealNeutrinoSolutions_ = mez.IsComplex() ? 0 : 2;
00243   const math::XYZTLorentzVector p4( ptr->px(), ptr->py(), pz, sqrt(ptr->px()*ptr->px() + ptr->py()*ptr->py() + pz*pz) );
00244   neutrino_ = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), p4, ptr->vertex() );
00245 }
00246 
00248 void
00249 TtSemiLepHypothesis::buildHypo(const edm::Handle<edm::View<reco::RecoCandidate> >& leps,
00250                                const edm::Handle<std::vector<pat::MET> >& mets, 
00251                                const edm::Handle<std::vector<pat::Jet> >& jets, 
00252                                std::vector<int>& match)
00253 {
00254   // -----------------------------------------------------
00255   // add jets
00256   // -----------------------------------------------------
00257   for(unsigned idx=0; idx<match.size(); ++idx){
00258     if( isValid(match[idx], jets) ){
00259       switch(idx){
00260       case TtSemiLepEvtPartons::LightQ:
00261         setCandidate(jets, match[idx], lightQ_,    jetCorrectionLevel("wQuarkMix")); break;
00262       case TtSemiLepEvtPartons::LightQBar:
00263         setCandidate(jets, match[idx], lightQBar_, jetCorrectionLevel("wQuarkMix")); break;
00264       case TtSemiLepEvtPartons::HadB:
00265         setCandidate(jets, match[idx], hadronicB_, jetCorrectionLevel("bQuark")); break;
00266       case TtSemiLepEvtPartons::LepB: 
00267         setCandidate(jets, match[idx], leptonicB_, jetCorrectionLevel("bQuark")); break;
00268       }
00269     }
00270   }
00271 
00272   // -----------------------------------------------------
00273   // add lepton
00274   // -----------------------------------------------------
00275   if( leps->empty() )
00276     return;
00277   setCandidate(leps, 0, lepton_);
00278   match.push_back( 0 );
00279   
00280   // -----------------------------------------------------
00281   // add neutrino
00282   // -----------------------------------------------------
00283   if( mets->empty() )
00284     return;
00285   if(neutrinoSolutionType_ == -1)
00286     setCandidate(mets, 0, neutrino_);
00287   else
00288     setNeutrino(mets, leps, 0, neutrinoSolutionType_);
00289 }
00290