CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepJetCombWMassMaxSumPt.cc

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