CMS 3D CMS Logo

TtSemiLepHypMaxSumPtWMass.cc

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

Generated on Tue Jun 9 17:48:13 2009 for CMSSW by  doxygen 1.5.4