CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepJetCombMaxSumPtWMass.cc

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