00001 #include "TopQuarkAnalysis/TopJetCombination/plugins/TtFullLepHypGenMatch.h"
00002 #include "AnalysisDataFormats/TopObjects/interface/TtFullLepEvtPartons.h"
00003 #include "DataFormats/Math/interface/deltaR.h"
00004
00005 TtFullLepHypGenMatch::TtFullLepHypGenMatch(const edm::ParameterSet& cfg):
00006 TtFullLepHypothesis( cfg )
00007 {
00008 }
00009
00010 TtFullLepHypGenMatch::~TtFullLepHypGenMatch() { }
00011
00012 void
00013 TtFullLepHypGenMatch::buildHypo(edm::Event& evt,
00014 const edm::Handle<std::vector<pat::Electron > >& elecs,
00015 const edm::Handle<std::vector<pat::Muon> >& mus,
00016 const edm::Handle<std::vector<pat::Jet> >& jets,
00017 const edm::Handle<std::vector<pat::MET> >& mets,
00018 std::vector<int>& match,
00019 const unsigned int iComb)
00020 {
00021
00022
00023
00024 for(unsigned idx=0; idx<match.size(); ++idx){
00025 if( isValid(match[idx], jets) ){
00026 switch(idx){
00027 case TtFullLepEvtPartons::B:
00028 setCandidate(jets, match[idx], b_ , jetCorrectionLevel_); break;
00029 case TtFullLepEvtPartons::BBar:
00030 setCandidate(jets, match[idx], bBar_, jetCorrectionLevel_); break;
00031 }
00032 }
00033 }
00034
00035
00036
00037
00038
00039 edm::Handle<TtGenEvent> genEvt;
00040 evt.getByLabel("genEvt", genEvt);
00041
00042
00043 if( !genEvt->isFullLeptonic() || !genEvt->lepton() || !genEvt->leptonBar() ){
00044 match.push_back( -1 );
00045 match.push_back( -1 );
00046 match.push_back( -1 );
00047 match.push_back( -1 );
00048 }
00049 else if(genEvt->isFullLeptonic(WDecay::kElec, WDecay::kElec) && elecs->size()>=2){
00050
00051 int iLepBar = findMatchingLepton(genEvt->leptonBar(), elecs);
00052 setCandidate(elecs, iLepBar, leptonBar_);
00053 match.push_back( iLepBar );
00054 int iLep = findMatchingLepton(genEvt->lepton(), elecs);
00055 setCandidate(elecs, iLep, lepton_);
00056 match.push_back( iLep );
00057
00058
00059 match.push_back( -1 );
00060 match.push_back( -1 );
00061 }
00062 else if(genEvt->isFullLeptonic(WDecay::kElec, WDecay::kMuon) && !elecs->empty() && !mus->empty() ){
00063 if(genEvt->leptonBar()->isElectron()){
00064
00065 int iLepBar = findMatchingLepton(genEvt->leptonBar(), elecs);
00066 setCandidate(elecs, iLepBar, leptonBar_);
00067 match.push_back( iLepBar );
00068
00069 match.push_back( -1 );
00070 match.push_back( -1 );
00071
00072 int iLep = findMatchingLepton(genEvt->lepton(), mus);
00073 setCandidate(mus, iLep, lepton_);
00074 match.push_back( iLep );
00075 }
00076 else{
00077
00078 match.push_back( -1 );
00079
00080 int iLepBar = findMatchingLepton(genEvt->leptonBar(), mus);
00081 setCandidate(mus, iLepBar, leptonBar_);
00082 match.push_back( iLepBar );
00083
00084 int iLep = findMatchingLepton(genEvt->lepton(), elecs);
00085 setCandidate(elecs, iLep, lepton_);
00086 match.push_back( iLep );
00087
00088 match.push_back( -1 );
00089 }
00090 }
00091 else if(genEvt->isFullLeptonic(WDecay::kMuon, WDecay::kMuon) && mus->size()>=2 ){
00092
00093 match.push_back( -1 );
00094 match.push_back( -1 );
00095
00096
00097 int iLepBar = findMatchingLepton(genEvt->leptonBar(), mus);
00098 setCandidate(mus, iLepBar, leptonBar_);
00099 match.push_back( iLepBar );
00100 int iLep = findMatchingLepton(genEvt->lepton(), mus);
00101 setCandidate(mus, iLep, lepton_);
00102 match.push_back( iLep );
00103 }
00104 else{
00105 match.push_back( -1 );
00106 match.push_back( -1 );
00107 match.push_back( -1 );
00108 match.push_back( -1 );
00109 }
00110
00111
00112
00113
00114 if( !mets->empty() ){
00115
00116 buildMatchingNeutrinos(evt, mets);
00117 }
00118 }
00119
00120 template <typename O>
00121 int
00122 TtFullLepHypGenMatch::findMatchingLepton(const reco::GenParticle* genLep, const edm::Handle<std::vector<O> >& leps)
00123 {
00124 int idx=-1;
00125 double minDR = -1;
00126 for(unsigned i=0; i<leps->size(); ++i){
00127 double dR = deltaR(genLep->eta(), genLep->phi(), (*leps)[i].eta(), (*leps)[i].phi());
00128 if(minDR<0 || dR<minDR){
00129 minDR=dR;
00130 idx=i;
00131 }
00132 }
00133 return idx;
00134 }
00135
00136 void
00137 TtFullLepHypGenMatch::buildMatchingNeutrinos(edm::Event& evt, const edm::Handle<std::vector<pat::MET> >& mets)
00138 {
00139
00140 edm::Handle<TtGenEvent> genEvt;
00141 evt.getByLabel("genEvt", genEvt);
00142
00143 if( genEvt->isTtBar() && genEvt->isFullLeptonic() && genEvt->neutrino() && genEvt->neutrinoBar() ){
00144 double momXNu = genEvt->neutrino() ->px();
00145 double momYNu = genEvt->neutrino() ->py();
00146 double momXNuBar = genEvt->neutrinoBar()->px();
00147 double momYNuBar = genEvt->neutrinoBar()->py();
00148
00149 double momXMet = mets->at(0).px();
00150 double momYMet = mets->at(0).py();
00151
00152 double momXNeutrino = 0.5*(momXNu - momXNuBar + momXMet);
00153 double momYNeutrino = 0.5*(momYNu - momYNuBar + momYMet);
00154 double momXNeutrinoBar = momXMet - momXNeutrino;
00155 double momYNeutrinoBar = momYMet - momYNeutrino;
00156
00157 math::XYZTLorentzVector recNuFM(momXNeutrino,momYNeutrino,0,sqrt(momXNeutrino * momXNeutrino + momYNeutrino * momYNeutrino));
00158 recNu = new reco::LeafCandidate(0,recNuFM);
00159
00160 math::XYZTLorentzVector recNuBarFM(momXNeutrinoBar,momYNeutrinoBar,0,sqrt(momXNeutrinoBar * momXNeutrinoBar + momYNeutrinoBar * momYNeutrinoBar));
00161 recNuBar = new reco::LeafCandidate(0,recNuBarFM);
00162 }
00163 }