00001 #include "FWCore/Framework/interface/Event.h"
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 #include "FWCore/ServiceRegistry/interface/Service.h"
00004 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006
00007 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLeptonicEvent.h"
00008
00009 #include "TopQuarkAnalysis/Examples/plugins/HypothesisAnalyzer.h"
00010
00011 HypothesisAnalyzer::HypothesisAnalyzer(const edm::ParameterSet& cfg):
00012 semiLepEvt_ (cfg.getParameter<edm::InputTag>("semiLepEvent")),
00013 hypoClassKey_(cfg.getParameter<edm::InputTag>("hypoClassKey"))
00014 {
00015 }
00016
00017 void
00018 HypothesisAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup)
00019 {
00021
00023
00024 edm::Handle<TtSemiLeptonicEvent> semiLepEvt;
00025 event.getByLabel(semiLepEvt_, semiLepEvt);
00026
00027 edm::Handle<int> hypoClassKeyHandle;
00028 event.getByLabel(hypoClassKey_, hypoClassKeyHandle);
00029 TtSemiLeptonicEvent::HypoClassKey& hypoClassKey = (TtSemiLeptonicEvent::HypoClassKey&) *hypoClassKeyHandle;
00030
00032
00034
00035 if( !semiLepEvt->isHypoValid(hypoClassKey) ){
00036 edm::LogInfo("HypothesisAnalyzer") << "Hypothesis not valid for this event";
00037 return;
00038 }
00039
00041
00043
00044 const reco::Candidate* hadTop = semiLepEvt->hadronicDecayTop(hypoClassKey);
00045 const reco::Candidate* hadW = semiLepEvt->hadronicDecayW (hypoClassKey);
00046
00048
00050
00051 if(hadW) {
00052 hadWPt_ ->Fill( hadW->pt() );
00053 hadWEta_ ->Fill( hadW->eta() );
00054 hadWMass_->Fill( hadW->mass() );
00055 }
00056
00057 if(hadTop) {
00058 hadTopPt_ ->Fill( hadTop->pt() );
00059 hadTopEta_ ->Fill( hadTop->eta() );
00060 hadTopMass_->Fill( hadTop->mass() );
00061 }
00062
00064
00066
00067 const reco::Candidate* genHadTop = semiLepEvt->hadronicDecayTop();
00068 const reco::Candidate* genHadW = semiLepEvt->hadronicDecayW();
00069
00071
00073
00074 if(hadW && genHadW) {
00075 hadWPullPt_ ->Fill( (hadW->pt() - genHadW->pt()) / genHadW->pt() );
00076 hadWPullEta_ ->Fill( (hadW->eta() - genHadW->eta()) / genHadW->eta() );
00077 hadWPullMass_->Fill( (hadW->mass() - genHadW->mass()) / genHadW->mass() );
00078 }
00079
00080 if(hadTop && genHadTop) {
00081 hadTopPullPt_ ->Fill( (hadTop->pt() - genHadTop->pt()) / genHadTop->pt() );
00082 hadTopPullEta_ ->Fill( (hadTop->eta() - genHadTop->eta()) / genHadTop->eta() );
00083 hadTopPullMass_->Fill( (hadTop->mass() - genHadTop->mass()) / genHadTop->mass() );
00084 }
00085
00087
00089
00090 genMatchDr_->Fill(semiLepEvt->genMatchSumDR());
00091 mvaDisc_ ->Fill(semiLepEvt->mvaDisc());
00092
00093 if(hadTop && genHadTop) {
00094
00095 genMatchDrVsHadTopPullMass_->Fill((hadTop->mass() - genHadTop->mass()) / genHadTop->mass(), semiLepEvt->genMatchSumDR());
00096 mvaDiscVsHadTopPullMass_ ->Fill((hadTop->mass() - genHadTop->mass()) / genHadTop->mass(), semiLepEvt->mvaDisc());
00097
00098 }
00099
00100 }
00101
00102 void
00103 HypothesisAnalyzer::beginJob()
00104 {
00105 edm::Service<TFileService> fs;
00106 if( !fs ) throw edm::Exception( edm::errors::Configuration, "TFile Service is not registered in cfg file" );
00107
00109
00111
00112 hadWPt_ = fs->make<TH1F>("hadWPt" , "p_{t} (W_{had}) [GeV]", 25, 0., 500.);
00113 hadWEta_ = fs->make<TH1F>("hadWEta" , "#eta (W_{had})" , 20, -4., 4.);
00114 hadWMass_ = fs->make<TH1F>("hadWMass", "M (W_{had}) [GeV]" , 25, 0., 200.);
00115
00116 hadTopPt_ = fs->make<TH1F>("hadTopPt" , "p_{t} (t_{had}) [GeV]", 25, 0. , 500.);
00117 hadTopEta_ = fs->make<TH1F>("hadTopEta" , "#eta (t_{had})" , 20, -4., 4.);
00118 hadTopMass_ = fs->make<TH1F>("hadTopMass", "M (t_{had}) [GeV]" , 40, 0. , 400.);
00119
00120 hadWPullPt_ = fs->make<TH1F>("hadWPullPt" , "(p_{t,rec}-p_{t,gen})/p_{t,gen} (W_{had})" , 40, -1., 1.);
00121 hadWPullEta_ = fs->make<TH1F>("hadWPullEta" , "(#eta_{rec}-#eta_{gen})/#eta_{gen} (W_{had})", 40, -1., 1.);
00122 hadWPullMass_ = fs->make<TH1F>("hadWPullMass", "(M_{rec}-M_{gen})/M_{gen} (W_{had})" , 40, -1., 1.);
00123
00124 hadTopPullPt_ = fs->make<TH1F>("hadTopPullPt" , "(p_{t,rec}-p_{t,gen})/p_{t,gen} (t_{had})" , 40, -1., 1.);
00125 hadTopPullEta_ = fs->make<TH1F>("hadTopPullEta" , "(#eta_{rec}-#eta_{gen})/#eta_{gen} (t_{had})", 40, -1., 1.);
00126 hadTopPullMass_ = fs->make<TH1F>("hadTopPullMass", "(M_{rec}-M_{gen})/M_{gen} (t_{had})" , 40, -1., 1.);
00127
00128 genMatchDr_ = fs->make<TH1F>("genMatchDr", "GenMatch #Sigma #Delta R", 40, 0., 4.);
00129 mvaDisc_ = fs->make<TH1F>("mvaDisc" , "MVA discriminator" , 20, 0., 1.);
00130
00131 genMatchDrVsHadTopPullMass_ = fs->make<TH2F>("genMatchDrVsHadTopPullMass",
00132 "GenMatch #Sigma #Delta R vs. (M_{rec}-M_{gen})/M_{gen} (t_{had}))",
00133 40, -1., 1., 40, 0., 4.);
00134 mvaDiscVsHadTopPullMass_ = fs->make<TH2F>("mvaDiscVsHadTopPullMass",
00135 "MVA discriminator vs. (M_{rec}-M_{gen})/M_{gen} (t_{had}))",
00136 40, -1., 1., 20, 0., 1.);
00137 }
00138
00139 void
00140 HypothesisAnalyzer::endJob()
00141 {
00142 }