CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/TopQuarkAnalysis/TopEventSelection/src/TtDilepLRSignalSelObservables.cc

Go to the documentation of this file.
00001 #include "TopQuarkAnalysis/TopEventSelection/interface/TtDilepLRSignalSelObservables.h"
00002 #include "DataFormats/PatCandidates/interface/Jet.h"
00003 #include "DataFormats/Math/interface/deltaR.h"
00004 #include "DataFormats/Math/interface/deltaPhi.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
00007 
00008 /************** Definition of the functions of the class ***************/
00009 
00010 //Constructor
00011 TtDilepLRSignalSelObservables::TtDilepLRSignalSelObservables(){
00012 count1=0; count2=0; count3=0;
00013 count4=0; count5=0; count3=0;
00014 }
00015 
00016 // Destructor
00017 TtDilepLRSignalSelObservables::~TtDilepLRSignalSelObservables(){
00018 }
00019 
00020 std::vector< TtDilepLRSignalSelObservables::IntBoolPair >
00021 TtDilepLRSignalSelObservables::operator() (TtDilepEvtSolution &solution,
00022         const edm::Event & iEvent, bool matchOnly)
00023 {
00024   evtselectVarVal.clear();
00025   evtselectVarMatch.clear();
00026   
00027   // Check whether the objects are matched:
00028   bool matchB1 = false;
00029   bool matchB2 = false;
00030   bool matchB = false;
00031   bool matchBbar = false;
00032   bool matchLeptPos =  false;
00033   bool matchLeptNeg =  false;
00034 
00035   try {
00036     // std::cout <<std::endl;
00037     double dr, dr1, dr2;
00038 
00039     edm::Handle<TtGenEvent> genEvent;
00040     iEvent.getByLabel ("genEvt",genEvent);
00041 
00042     if (genEvent->isFullLeptonic()) {
00043       // Match the leptons, by type and deltaR
00044       dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptPos(), *(solution.getGenLepp()));
00045       matchLeptPos = (
00046         ( ((solution.getWpDecay()=="electron")&&(std::abs(solution.getGenLepp()->pdgId())==11))
00047        || ((solution.getWpDecay()=="muon")&&(std::abs(solution.getGenLepp()->pdgId())==13)) )
00048        && (dr < 0.1) );
00049 
00050       dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptNeg(), *(solution.getGenLepm()));
00051       matchLeptNeg = (
00052         ( ((solution.getWmDecay()=="electron")&&(std::abs(solution.getGenLepm()->pdgId())==11))
00053            || ((solution.getWmDecay()=="muon")&&(std::abs(solution.getGenLepm()->pdgId())==13)) )
00054         && (dr < 0.1) );
00055     }
00056 
00057     if (genEvent->isSemiLeptonic()) {
00058       int id = genEvent->singleLepton()->pdgId();
00059 
00060       if (id>0) {
00061         dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptNeg(), *(genEvent->singleLepton()));
00062         matchLeptNeg = (
00063           ( ((solution.getWmDecay()=="electron") && (id==11))
00064              || ((solution.getWmDecay()=="muon") && (id==13)) )
00065           && (dr < 0.1) );
00066       } else {
00067         dr = DeltaR<reco::Particle, reco::GenParticle>()(solution.getLeptPos(), *(genEvent->singleLepton()));
00068         matchLeptPos = (
00069           ( ((solution.getWpDecay()=="electron")&& (id==-11))
00070          || ((solution.getWpDecay()=="muon")    && (id==-13)) )
00071          && (dr < 0.1) );
00072       }
00073     }
00074     
00075     if (genEvent->isTtBar() && genEvent->numberOfBQuarks()>1) {
00076       if (solution.getJetB().partonFlavour()==5) ++count1;
00077       if (solution.getJetBbar().partonFlavour()==5) ++count1;
00078 
00079       dr1 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetB(), *(genEvent->b()));
00080       dr2 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetB(), *(genEvent->bBar()));
00081 
00082       matchB1= ( (dr1<0.4) || (dr2<0.4));
00083       matchB = ( (solution.getJetB().partonFlavour()==5) && (dr1<0.4) );
00084       if (matchB) ++count3;
00085       matchB = ( (dr1<0.4) );
00086       if (dr1<0.5) ++count2;
00087       if (dr1<0.4) ++count4;
00088       if (dr1<0.3) ++count5;
00089 
00090       dr1 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetBbar(), *(genEvent->b()));
00091       dr2 = DeltaR<pat::Jet, reco::GenParticle>()(solution.getCalJetBbar(), *(genEvent->bBar()));
00092 
00093       matchBbar = ( (solution.getJetBbar().partonFlavour()==5) && (dr2<0.4) );
00094       if (matchBbar) ++count3;
00095       matchBbar = ( (dr2<0.4) );
00096       matchB2 = ( (dr1<0.4) || (dr2<0.4));
00097       if (dr2<0.5) ++count2;
00098       if (dr2<0.4) ++count4;
00099       if (dr2<0.3) ++count5;
00100     }
00101     
00102   } catch (...){std::cout << "Exception\n";}
00103   
00104   edm::Handle<std::vector<pat::Jet> > jets;
00105   iEvent.getByLabel(jetSource_, jets);
00106   
00107   //  Lower / Higher of both jet angles
00108   
00109   double v1 = std::abs( solution.getJetB().p4().theta() - M_PI/2 );
00110   double v2 = std::abs( solution.getJetBbar().p4().theta() - M_PI/2 ) ;
00111   fillMinMax(v1, v2, 1, evtselectVarVal, matchB1, matchB2, evtselectVarMatch);
00112 
00113   //  Lower / Higher of both jet pT
00114 
00115   double pt1 = solution.getJetB().p4().pt();
00116   double pt2 = solution.getJetBbar().p4().pt();
00117   fillMinMax(pt1, pt2, 3, evtselectVarVal, matchB1, matchB2, evtselectVarMatch);
00118 
00119   //  Lower / Higher of both lepton pT
00120 
00121   pt1 = solution.getLeptPos().p4().pt();
00122   pt2 = solution.getLeptNeg().p4().pt();
00123   fillMinMax(pt1, pt2, 5, evtselectVarVal, matchLeptPos, matchLeptNeg, evtselectVarMatch);
00124 
00125   // delta theta btw the b-jets
00126 
00127   double deltaPhi = std::abs( delta(solution.getJetB().p4().phi(), 
00128                                     solution.getJetBbar().p4().phi()) );
00129 
00130   evtselectVarVal.push_back(IntDblPair(7, deltaPhi));
00131   evtselectVarMatch.push_back(IntBoolPair(7, matchB1&&matchB2));
00132 
00133   // delta phi btw the b-jets
00134 
00135   double deltaTheta = std::abs( delta (solution.getJetBbar().p4().theta(),
00136                                        solution.getJetB().p4().theta() ) );
00137 
00138   evtselectVarVal.push_back(IntDblPair(8, deltaTheta));
00139   evtselectVarMatch.push_back(IntBoolPair(8, matchB1&&matchB2));
00140 
00141   //  Lower / Higher of phi difference between the b and associated lepton
00142 
00143   double deltaPhi1 = std::abs ( delta( solution.getJetB().p4().phi(),
00144                                        solution.getLeptPos().p4().phi() ) );
00145   double deltaPhi2 = std::abs ( delta( solution.getJetBbar().p4().phi(),
00146                                        solution.getLeptNeg().p4().phi() ) );
00147 
00148   fillMinMax(deltaPhi1, deltaPhi2, 9, evtselectVarVal, 
00149         matchB&&matchLeptPos, matchBbar&&matchLeptNeg, evtselectVarMatch);
00150 
00151   //  Lower / Higher of theta difference between the b and associated lepton
00152 
00153   double deltaTheta1 = std::abs( solution.getJetB().p4().theta() -
00154                                  solution.getLeptPos().p4().theta() );
00155   double deltaTheta2 = std::abs( solution.getJetBbar().p4().theta() -
00156                                  solution.getLeptNeg().p4().theta() );
00157   fillMinMax(deltaTheta1, deltaTheta2, 11, evtselectVarVal, 
00158         matchB&&matchLeptPos, matchBbar&&matchLeptNeg, evtselectVarMatch);
00159 
00160   // Invariant Mass of lepton pair
00161 
00162   math::XYZTLorentzVector pp = solution.getLeptPos().p4() + solution.getLeptNeg().p4();
00163   double mass = pp.mass();
00164   evtselectVarVal.push_back(IntDblPair(13, mass));
00165   evtselectVarMatch.push_back(IntBoolPair(13, matchLeptNeg&&matchLeptPos));
00166 
00167   evtselectVarVal.push_back(IntDblPair(13, mass));
00168   evtselectVarMatch.push_back(IntBoolPair(13, matchLeptNeg&&matchLeptPos));
00169 
00170   std::vector <pat::Jet> jet3;
00171   for (int i=0;i<(int)jets->size();++i) {
00172 if  ( ((*jets)[i].et()<solution.getJetB().et()) && ((*jets)[i].et()<solution.getJetBbar().et())) {jet3.push_back((*jets)[i]);
00173 }}
00174   double jet1Ratio = 0., jet2Ratio = 0.;  
00175   if (jet3.size()>0) { 
00176     jet1Ratio = jet3[0].et()/solution.getJetB().et();
00177     jet2Ratio = jet3[0].et()/solution.getJetBbar().et();
00178   }
00179   fillMinMax(jet1Ratio, jet2Ratio, 14, evtselectVarVal, 
00180         matchB1, matchB2, evtselectVarMatch);
00181 
00182   evtselectVarVal.push_back(IntDblPair(16, jets->size()));
00183   evtselectVarMatch.push_back(IntBoolPair(16, matchB&&matchBbar));
00184 
00185 
00186   if (!matchOnly) solution.setLRSignalEvtObservables(evtselectVarVal);
00187   return evtselectVarMatch;
00188 }
00189 
00190 void TtDilepLRSignalSelObservables::fillMinMax
00191         (double v1, double v2, int obsNbr, std::vector< IntDblPair > & varList,
00192          bool match1, bool match2, std::vector< IntBoolPair > & matchList)
00193 {
00194   if (v1<v2) {
00195     varList.push_back(IntDblPair(obsNbr, v1));
00196     varList.push_back(IntDblPair(obsNbr+1, v2));
00197     matchList.push_back(IntBoolPair(obsNbr, match1));
00198     matchList.push_back(IntBoolPair(obsNbr+1, match2));
00199     
00200   } else {
00201     varList.push_back(IntDblPair(obsNbr, v2));
00202     varList.push_back(IntDblPair(obsNbr+1, v1));
00203     matchList.push_back(IntBoolPair(obsNbr, match2));
00204     matchList.push_back(IntBoolPair(obsNbr+1, match1));
00205   }
00206 }
00207 
00208 double TtDilepLRSignalSelObservables::delta(double phi1, double phi2)
00209 {
00210   double deltaPhi = phi1 - phi2;
00211   while (deltaPhi > M_PI) deltaPhi -= 2*M_PI;
00212   while (deltaPhi <= -M_PI) deltaPhi += 2*M_PI;
00213   return deltaPhi;
00214 }