00001 #include "DataFormats/Candidate/interface/Candidate.h" 00002 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00003 #include "AnalysisDataFormats/TopObjects/interface/TtEvent.h" 00004 #include "TopQuarkAnalysis/Examples/plugins/HypothesisAnalyzer.h" 00005 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLeptonicEvent.h" 00006 00007 00008 HypothesisAnalyzer::HypothesisAnalyzer(const edm::ParameterSet& cfg): 00009 semiLepEvt_ (cfg.getParameter<edm::InputTag>("semiLepEvent")), 00010 hypoClassKey_(cfg.getParameter<edm::InputTag>("hypoClassKey")) 00011 { 00012 } 00013 00014 void 00015 HypothesisAnalyzer::analyze(const edm::Event& evt, const edm::EventSetup& setup) 00016 { 00017 edm::Handle<TtSemiLeptonicEvent> semiLepEvt; 00018 evt.getByLabel(semiLepEvt_, semiLepEvt); 00019 00020 edm::Handle<int> hypoClassKeyHandle; 00021 evt.getByLabel(hypoClassKey_, hypoClassKeyHandle); 00022 TtSemiLeptonicEvent::HypoClassKey& hypoClassKey = (TtSemiLeptonicEvent::HypoClassKey&) *hypoClassKeyHandle; 00023 00024 if( !semiLepEvt->isHypoAvailable(hypoClassKey) ){ 00025 edm::LogInfo ( "NonValidHyp" ) << "Hypothesis not available for this event"; 00026 return; 00027 } 00028 if( !semiLepEvt->isHypoValid(hypoClassKey) ){ 00029 edm::LogInfo ( "NonValidHyp" ) << "Hypothesis not valid for this event"; 00030 return; 00031 } 00032 00033 const reco::Candidate* hadTop = semiLepEvt->hadronicTop(hypoClassKey); 00034 const reco::Candidate* hadW = semiLepEvt->hadronicW (hypoClassKey); 00035 const reco::Candidate* lepTop = semiLepEvt->leptonicTop(hypoClassKey); 00036 const reco::Candidate* lepW = semiLepEvt->leptonicW (hypoClassKey); 00037 00038 if(hadTop && hadW && lepTop && lepW){ 00039 hadWPt_ ->Fill( hadW->pt() ); 00040 hadWMass_ ->Fill( hadW->mass() ); 00041 hadTopPt_ ->Fill( hadTop->pt() ); 00042 hadTopMass_->Fill( hadTop->mass()); 00043 00044 lepWPt_ ->Fill( lepW->pt() ); 00045 lepWMass_ ->Fill( lepW->mass() ); 00046 lepTopPt_ ->Fill( lepTop->pt() ); 00047 lepTopMass_->Fill( lepTop->mass()); 00048 } 00049 } 00050 00051 void 00052 HypothesisAnalyzer::beginJob(const edm::EventSetup&) 00053 { 00054 edm::Service<TFileService> fs; 00055 if( !fs ) throw edm::Exception( edm::errors::Configuration, "TFile Service is not registered in cfg file" ); 00056 00057 hadWPt_ = fs->make<TH1F>("hadWPt", "p_{t} (W_{had}) [GeV]", 100, 0., 500.); 00058 hadWMass_ = fs->make<TH1F>("hadWMass", "M (W_{had}) [GeV]" , 50, 0. , 150.); 00059 hadTopPt_ = fs->make<TH1F>("hadTopPt", "p_{t} (t_{had}) [GeV]", 100, 0. , 500.); 00060 hadTopMass_ = fs->make<TH1F>("hadTopMass", "M (t_{had}) [GeV]", 50, 50. , 250.); 00061 00062 lepWPt_ = fs->make<TH1F>("lepWPt", "p_{t} (W_{lep}) [GeV]", 100, 0., 500.); 00063 lepWMass_ = fs->make<TH1F>("lepWMass", "M (W_{lep}) [GeV]" , 50, 0. , 150.); 00064 lepTopPt_ = fs->make<TH1F>("lepTopPt", "p_{t} (t_{lep}) [GeV]", 100, 0. , 500.); 00065 lepTopMass_ = fs->make<TH1F>("lepTopMass", "M (t_{lep}) [GeV]", 50, 50. , 250.); 00066 } 00067 00068 void 00069 HypothesisAnalyzer::endJob() 00070 { 00071 }