CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepHypWMassDeltaTopMass.cc

Go to the documentation of this file.
00001 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00002 #include "TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepHypWMassDeltaTopMass.h"
00003 
00004 TtSemiLepHypWMassDeltaTopMass::TtSemiLepHypWMassDeltaTopMass(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 TtSemiLepHypWMassDeltaTopMass::~TtSemiLepHypWMassDeltaTopMass() { }
00021 
00022 void
00023 TtSemiLepHypWMassDeltaTopMass::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   // add lepton
00054   // -----------------------------------------------------
00055   setCandidate(leps, 0, lepton_);
00056   match[TtSemiLepEvtPartons::Lepton] = 0;
00057   
00058   // -----------------------------------------------------
00059   // add neutrino
00060   // -----------------------------------------------------
00061   if(neutrinoSolutionType_ == -1)
00062     setCandidate(mets, 0, neutrino_);
00063   else
00064     setNeutrino(mets, leps, 0, neutrinoSolutionType_);
00065 
00066   // -----------------------------------------------------
00067   // associate those jets that get closest to the W mass
00068   // with their invariant mass to the hadronic W boson
00069   // -----------------------------------------------------
00070   double wDist =-1.;
00071   std::vector<int> closestToWMassIndices;
00072   closestToWMassIndices.push_back(-1);
00073   closestToWMassIndices.push_back(-1);
00074   for(int idx=0; idx<maxNJets; ++idx){
00075     if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
00076     for(int jdx=(idx+1); jdx<maxNJets; ++jdx){
00077       if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
00078       reco::Particle::LorentzVector sum = 
00079         (*jets)[idx].p4()+
00080         (*jets)[jdx].p4();
00081       if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
00082         wDist=fabs(sum.mass()-wMass_);
00083         closestToWMassIndices.clear();
00084         closestToWMassIndices.push_back(idx);
00085         closestToWMassIndices.push_back(jdx);
00086       }
00087     }
00088   }
00089 
00090   // -----------------------------------------------------
00091   // associate those jets to the hadronic and the leptonic
00092   // b quark that minimize the difference between
00093   // hadronic and leptonic top mass
00094   // -----------------------------------------------------
00095   double deltaTop=-1.;
00096   int hadB=-1;
00097   int lepB=-1;
00098   if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
00099     reco::Particle::LorentzVector lepW = neutrino_->p4() + lepton_->p4();
00100     reco::Particle::LorentzVector hadW =
00101       (*jets)[closestToWMassIndices[0]].p4()+
00102       (*jets)[closestToWMassIndices[1]].p4();
00103     // find hadronic b candidate
00104     for(int idx=0; idx<maxNJets; ++idx){
00105       if(useBTagging_ && !isBJet[idx]) continue;
00106       // make sure it's not used up already from the hadronic W
00107       if( idx!=closestToWMassIndices[0] && idx!=closestToWMassIndices[1] ){
00108         reco::Particle::LorentzVector hadTop = hadW + (*jets)[idx].p4();
00109         // find leptonic b candidate
00110         for(int jdx=0; jdx<maxNJets; ++jdx){
00111           if(useBTagging_ && !isBJet[jdx]) continue;
00112           // make sure it's not used up already from the hadronic branch
00113           if( jdx!=closestToWMassIndices[0] && jdx!=closestToWMassIndices[1] && jdx!=idx ){
00114             reco::Particle::LorentzVector lepTop = lepW + (*jets)[jdx].p4();
00115             if( deltaTop<0. || deltaTop>fabs(hadTop.mass()-lepTop.mass()) ){
00116               deltaTop=fabs(hadTop.mass()-lepTop.mass());
00117               hadB=idx;
00118               lepB=jdx;
00119             }
00120           }
00121         }
00122       }
00123     }
00124   }
00125 
00126   // -----------------------------------------------------
00127   // add jets
00128   // -----------------------------------------------------
00129   if( isValid(closestToWMassIndices[0], jets) ){
00130     setCandidate(jets, closestToWMassIndices[0], lightQ_, jetCorrectionLevel("wQuarkMix"));
00131     match[TtSemiLepEvtPartons::LightQ] = closestToWMassIndices[0];
00132   }
00133 
00134   if( isValid(closestToWMassIndices[1], jets) ){
00135     setCandidate(jets, closestToWMassIndices[1], lightQBar_, jetCorrectionLevel("wQuarkMix"));
00136     match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
00137   }
00138   
00139   if( isValid(hadB, jets) ){
00140     setCandidate(jets, hadB, hadronicB_, jetCorrectionLevel("bQuark"));
00141     match[TtSemiLepEvtPartons::HadB] = hadB;
00142   }
00143 
00144   if( isValid(lepB, jets) ){
00145     setCandidate(jets, lepB, leptonicB_, jetCorrectionLevel("bQuark"));
00146     match[TtSemiLepEvtPartons::LepB] = lepB;
00147   }
00148 
00149 }