CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepHypWMassMaxSumPt.cc

Go to the documentation of this file.
00001 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00002 #include "TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepHypWMassMaxSumPt.h"
00003 
00004 TtSemiLepHypWMassMaxSumPt::TtSemiLepHypWMassMaxSumPt(const edm::ParameterSet& cfg):
00005   TtSemiLepHypothesis( cfg ),  
00006   maxNJets_            (cfg.getParameter<int>        ("maxNJets"            )),
00007   wMass_               (cfg.getParameter<double>     ("wMass"               )),
00008   useBTagging_         (cfg.getParameter<bool>       ("useBTagging"         )),
00009   bTagAlgorithm_       (cfg.getParameter<std::string>("bTagAlgorithm"       )),
00010   minBDiscBJets_       (cfg.getParameter<double>     ("minBDiscBJets"       )),
00011   maxBDiscLightJets_   (cfg.getParameter<double>     ("maxBDiscLightJets"   )),
00012   neutrinoSolutionType_(cfg.getParameter<int>        ("neutrinoSolutionType"))
00013 {
00014   if(maxNJets_<4 && maxNJets_!=-1)
00015     throw cms::Exception("WrongConfig") 
00016       << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
00017       << "It has to be larger than 4 or can be set to -1 to take all jets.";
00018 }
00019 
00020 TtSemiLepHypWMassMaxSumPt::~TtSemiLepHypWMassMaxSumPt() { }
00021 
00022 void
00023 TtSemiLepHypWMassMaxSumPt::buildHypo(edm::Event& evt,
00024                                      const edm::Handle<edm::View<reco::RecoCandidate> >& leps,
00025                                      const edm::Handle<std::vector<pat::MET> >& mets, 
00026                                      const edm::Handle<std::vector<pat::Jet> >& jets, 
00027                                      std::vector<int>& match, const unsigned int iComb)
00028 {
00029   if(leps->empty() || mets->empty() || jets->size()<4){
00030     // create empty hypothesis
00031     return;
00032   }
00033 
00034   int maxNJets = maxNJets_;
00035   if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
00036 
00037   std::vector<bool> isBJet;
00038   std::vector<bool> isLJet;
00039   int cntBJets = 0;
00040   if(useBTagging_) {
00041     for(int idx=0; idx<maxNJets; ++idx) {
00042       isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_    ) );
00043       isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
00044       if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_    )cntBJets++;
00045     }
00046   }
00047 
00048   match.clear();
00049   for(unsigned int i=0; i<5; ++i)
00050     match.push_back(-1);
00051 
00052   // -----------------------------------------------------
00053   // associate those jets that get closest to the W mass
00054   // with their invariant mass to the hadronic W boson
00055   // -----------------------------------------------------
00056   double wDist =-1.;
00057   std::vector<int> closestToWMassIndices;
00058   closestToWMassIndices.push_back(-1);
00059   closestToWMassIndices.push_back(-1);
00060   for(int idx=0; idx<maxNJets; ++idx){
00061     if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
00062     for(int jdx=(idx+1); jdx<maxNJets; ++jdx){
00063       if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
00064       reco::Particle::LorentzVector sum = 
00065         (*jets)[idx].p4()+
00066         (*jets)[jdx].p4();
00067       if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
00068         wDist=fabs(sum.mass()-wMass_);
00069         closestToWMassIndices.clear();
00070         closestToWMassIndices.push_back(idx);
00071         closestToWMassIndices.push_back(jdx);
00072       }
00073     }
00074   }
00075 
00076   // -----------------------------------------------------
00077   // associate those jets with maximum pt of the vectorial 
00078   // sum to the hadronic decay chain
00079   // -----------------------------------------------------
00080   double maxPt=-1.;
00081   int hadB=-1;
00082   if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
00083     for(int idx=0; idx<maxNJets; ++idx){
00084       if(useBTagging_ && !isBJet[idx]) continue;
00085       // make sure it's not used up already from the hadronic W
00086       if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] ){
00087         reco::Particle::LorentzVector sum = 
00088           (*jets)[closestToWMassIndices[0]].p4()+
00089           (*jets)[closestToWMassIndices[1]].p4()+
00090           (*jets)[idx].p4();
00091         if( maxPt<0. || maxPt<sum.pt() ){
00092           maxPt=sum.pt();
00093           hadB=idx;
00094         }
00095       }
00096     }
00097   }
00098   
00099   // -----------------------------------------------------
00100   // associate the remaining jet with maximum pt of the   
00101   // vectorial sum with the leading lepton with the 
00102   // leptonic b quark
00103   // -----------------------------------------------------
00104   maxPt=-1.;
00105   int lepB=-1;
00106   for(int idx=0; idx<maxNJets; ++idx){
00107     if(useBTagging_ && !isBJet[idx]) continue;
00108     // make sure it's not used up already from the hadronic decay chain
00109     if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] && idx!=hadB) {
00110       reco::Particle::LorentzVector sum = 
00111         (*jets)[idx].p4()+(*leps)[ 0 ].p4();
00112       if( maxPt<0. || maxPt<sum.pt() ){
00113         maxPt=sum.pt();
00114         lepB=idx;
00115       }
00116     }
00117   }
00118 
00119   // -----------------------------------------------------
00120   // add jets
00121   // -----------------------------------------------------
00122   if( isValid(closestToWMassIndices[0], jets) ){
00123     setCandidate(jets, closestToWMassIndices[0], lightQ_, jetCorrectionLevel("wQuarkMix"));
00124     match[TtSemiLepEvtPartons::LightQ] = closestToWMassIndices[0];
00125   }
00126 
00127   if( isValid(closestToWMassIndices[1], jets) ){
00128     setCandidate(jets, closestToWMassIndices[1], lightQBar_, jetCorrectionLevel("wQuarkMix"));
00129     match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
00130   }
00131   
00132   if( isValid(hadB, jets) ){
00133     setCandidate(jets, hadB, hadronicB_, jetCorrectionLevel("bQuark"));
00134     match[TtSemiLepEvtPartons::HadB] = hadB;
00135   }
00136 
00137   if( isValid(lepB, jets) ){
00138     setCandidate(jets, lepB, leptonicB_, jetCorrectionLevel("bQuark"));
00139     match[TtSemiLepEvtPartons::LepB] = lepB;
00140   }
00141 
00142   // -----------------------------------------------------
00143   // add lepton
00144   // -----------------------------------------------------
00145   setCandidate(leps, 0, lepton_);
00146   match[TtSemiLepEvtPartons::Lepton] = 0;
00147   
00148   // -----------------------------------------------------
00149   // add neutrino
00150   // -----------------------------------------------------
00151   if(neutrinoSolutionType_ == -1)
00152     setCandidate(mets, 0, neutrino_);
00153   else
00154     setNeutrino(mets, leps, 0, neutrinoSolutionType_);
00155 }