#include <HypothesisAnalyzer.h>
Public Member Functions | |
HypothesisAnalyzer (const edm::ParameterSet &) | |
~HypothesisAnalyzer () | |
Private Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
virtual void | endJob () |
Private Attributes | |
TH1F * | genMatchDr_ |
TH2F * | genMatchDrVsHadTopPullMass_ |
TH1F * | hadTopEta_ |
TH1F * | hadTopMass_ |
TH1F * | hadTopPt_ |
TH1F * | hadTopPullEta_ |
TH1F * | hadTopPullMass_ |
TH1F * | hadTopPullPt_ |
TH1F * | hadWEta_ |
TH1F * | hadWMass_ |
TH1F * | hadWPt_ |
TH1F * | hadWPullEta_ |
TH1F * | hadWPullMass_ |
TH1F * | hadWPullPt_ |
const std::string | hypoClassKey_ |
TH1F * | kinFitProb_ |
TH2F * | kinFitProbVsHadTopPullMass_ |
TH1F * | lepTopEta_ |
TH1F * | lepTopMass_ |
TH1F * | lepTopPt_ |
TH1F * | lepTopPullEta_ |
TH1F * | lepTopPullMass_ |
TH1F * | lepTopPullPt_ |
TH1F * | lepWEta_ |
TH1F * | lepWMass_ |
TH1F * | lepWPt_ |
TH1F * | lepWPullEta_ |
TH1F * | lepWPullMass_ |
TH1F * | lepWPullPt_ |
TH1F * | neutrinoEta_ |
TH1F * | neutrinoPullEta_ |
const edm::InputTag | semiLepEvt_ |
TH1F * | topPairMass_ |
TH1F * | topPairPullMass_ |
Definition at line 9 of file HypothesisAnalyzer.h.
HypothesisAnalyzer::HypothesisAnalyzer | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 11 of file HypothesisAnalyzer.cc.
: semiLepEvt_ (cfg.getParameter<edm::InputTag>("semiLepEvent")), hypoClassKey_(cfg.getParameter<std::string>("hypoClassKey")) { }
HypothesisAnalyzer::~HypothesisAnalyzer | ( | ) | [inline] |
Definition at line 14 of file HypothesisAnalyzer.h.
{};
void HypothesisAnalyzer::analyze | ( | const edm::Event & | event, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
Implements edm::EDAnalyzer.
Definition at line 18 of file HypothesisAnalyzer.cc.
References reco::Candidate::eta(), genMatchDr_, genMatchDrVsHadTopPullMass_, hadTopEta_, hadTopMass_, hadTopPt_, hadTopPullEta_, hadTopPullMass_, hadTopPullPt_, hadWEta_, hadWMass_, hadWPt_, hadWPullEta_, hadWPullMass_, hadWPullPt_, hypoClassKey_, kinFitProb_, kinFitProbVsHadTopPullMass_, lepTopEta_, lepTopMass_, lepTopPt_, lepTopPullEta_, lepTopPullMass_, lepTopPullPt_, lepWEta_, lepWMass_, lepWPt_, lepWPullEta_, lepWPullMass_, lepWPullPt_, reco::Candidate::mass(), neutrinoEta_, neutrinoPullEta_, reco::Candidate::pt(), semiLepEvt_, topPairMass_, and topPairPullMass_.
{ // get a handle for the TtSemiLeptonicEvent and a key to the hypothesis edm::Handle<TtSemiLeptonicEvent> semiLepEvt; event.getByLabel(semiLepEvt_, semiLepEvt); // check if hypothesis is available and valid in this event if( !semiLepEvt->isHypoValid(hypoClassKey_) ){ edm::LogInfo("HypothesisAnalyzer") << "Hypothesis " << hypoClassKey_ << " not valid for this event"; return; } // get reconstructed top quarks, W bosons, the top pair and the neutrino from the hypothesis const reco::Candidate* topPair = semiLepEvt->topPair(hypoClassKey_); const reco::Candidate* lepTop = semiLepEvt->leptonicDecayTop(hypoClassKey_); const reco::Candidate* lepW = semiLepEvt->leptonicDecayW(hypoClassKey_); const reco::Candidate* hadTop = semiLepEvt->hadronicDecayTop(hypoClassKey_); const reco::Candidate* hadW = semiLepEvt->hadronicDecayW(hypoClassKey_); const reco::Candidate* neutrino = semiLepEvt->singleNeutrino(hypoClassKey_); // fill simple histograms with kinematic variables of the reconstructed particles if(topPair) topPairMass_->Fill( topPair->mass() ); if(hadW) { hadWPt_ ->Fill( hadW->pt() ); hadWEta_ ->Fill( hadW->eta() ); hadWMass_->Fill( hadW->mass() ); } if(hadTop) { hadTopPt_ ->Fill( hadTop->pt() ); hadTopEta_ ->Fill( hadTop->eta() ); hadTopMass_->Fill( hadTop->mass() ); } if(lepW) { lepWPt_ ->Fill( lepW->pt() ); lepWEta_ ->Fill( lepW->eta() ); lepWMass_->Fill( lepW->mass() ); } if(lepTop) { lepTopPt_ ->Fill( lepTop->pt() ); lepTopEta_ ->Fill( lepTop->eta() ); lepTopMass_->Fill( lepTop->mass() ); } if(neutrino) neutrinoEta_->Fill( neutrino->eta() ); // get corresponding genParticles const math::XYZTLorentzVector* genTopPair = semiLepEvt->topPair(); const reco::Candidate* genHadTop = semiLepEvt->hadronicDecayTop(); const reco::Candidate* genHadW = semiLepEvt->hadronicDecayW(); const reco::Candidate* genLepTop = semiLepEvt->leptonicDecayTop(); const reco::Candidate* genLepW = semiLepEvt->leptonicDecayW(); const reco::Candidate* genNeutrino = semiLepEvt->singleNeutrino(); // fill pull histograms of kinematic variables with respect to the generated particles if(topPair && genTopPair) topPairPullMass_->Fill( (topPair->mass()-genTopPair->mass())/ genTopPair->mass() ); if(hadW && genHadW) { hadWPullPt_ ->Fill( (hadW->pt() - genHadW->pt()) / genHadW->pt() ); hadWPullEta_ ->Fill( (hadW->eta() - genHadW->eta()) / genHadW->eta() ); hadWPullMass_->Fill( (hadW->mass() - genHadW->mass()) / genHadW->mass() ); } if(hadTop && genHadTop) { hadTopPullPt_ ->Fill( (hadTop->pt() - genHadTop->pt()) / genHadTop->pt() ); hadTopPullEta_ ->Fill( (hadTop->eta() - genHadTop->eta()) / genHadTop->eta() ); hadTopPullMass_->Fill( (hadTop->mass() - genHadTop->mass()) / genHadTop->mass() ); } if(lepW && genLepW) { lepWPullPt_ ->Fill( (lepW->pt() - genLepW->pt()) / genLepW->pt() ); lepWPullEta_ ->Fill( (lepW->eta() - genLepW->eta()) / genLepW->eta() ); lepWPullMass_->Fill( (lepW->mass() - genLepW->mass()) / genLepW->mass() ); } if(lepTop && genLepTop) { lepTopPullPt_ ->Fill( (lepTop->pt() - genLepTop->pt()) / genLepTop->pt() ); lepTopPullEta_ ->Fill( (lepTop->eta() - genLepTop->eta()) / genLepTop->eta() ); lepTopPullMass_->Fill( (lepTop->mass() - genLepTop->mass()) / genLepTop->mass() ); } if(neutrino && genNeutrino) neutrinoPullEta_->Fill( (neutrino->eta()-genNeutrino->eta()) / genNeutrino->eta() ); // fill histograms with variables describing the quality of the hypotheses genMatchDr_->Fill(semiLepEvt->genMatchSumDR()); kinFitProb_->Fill(semiLepEvt->fitProb()); if(hadTop && genHadTop) { genMatchDrVsHadTopPullMass_->Fill((hadTop->mass() - genHadTop->mass()) / genHadTop->mass(), semiLepEvt->genMatchSumDR()); kinFitProbVsHadTopPullMass_->Fill((hadTop->mass() - genHadTop->mass()) / genHadTop->mass(), semiLepEvt->fitProb()); } }
void HypothesisAnalyzer::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 133 of file HypothesisAnalyzer.cc.
References edm::errors::Configuration, Exception, genMatchDr_, genMatchDrVsHadTopPullMass_, hadTopEta_, hadTopMass_, hadTopPt_, hadTopPullEta_, hadTopPullMass_, hadTopPullPt_, hadWEta_, hadWMass_, hadWPt_, hadWPullEta_, hadWPullMass_, hadWPullPt_, kinFitProb_, kinFitProbVsHadTopPullMass_, lepTopEta_, lepTopMass_, lepTopPt_, lepTopPullEta_, lepTopPullMass_, lepTopPullPt_, lepWEta_, lepWMass_, lepWPt_, lepWPullEta_, lepWPullMass_, lepWPullPt_, neutrinoEta_, neutrinoPullEta_, topPairMass_, and topPairPullMass_.
{ edm::Service<TFileService> fs; if( !fs ) throw edm::Exception( edm::errors::Configuration, "TFile Service is not registered in cfg file" ); // book histograms neutrinoEta_ = fs->make<TH1F>("neutrinoEta", "#eta (neutrino)", 21, -4., 4.); neutrinoPullEta_ = fs->make<TH1F>("neutrinoPullEta", "(#eta_{rec}-#eta_{gen})/#eta_{gen} (neutrino)", 40, -1., 1.); hadWPt_ = fs->make<TH1F>("hadWPt" , "p_{T} (W_{had}) [GeV]", 25, 0., 500.); hadWEta_ = fs->make<TH1F>("hadWEta" , "#eta (W_{had})" , 21, -4., 4.); hadWMass_ = fs->make<TH1F>("hadWMass", "M (W_{had}) [GeV]" , 25, 0., 200.); hadTopPt_ = fs->make<TH1F>("hadTopPt" , "p_{T} (t_{had}) [GeV]", 25, 0. , 500.); hadTopEta_ = fs->make<TH1F>("hadTopEta" , "#eta (t_{had})" , 21, -4., 4.); hadTopMass_ = fs->make<TH1F>("hadTopMass", "M (t_{had}) [GeV]" , 40, 0. , 400.); lepWPt_ = fs->make<TH1F>("lepWPt" , "p_{t} (W_{lep}) [GeV]", 25, 0., 500.); lepWEta_ = fs->make<TH1F>("lepWEta" , "#eta (W_{lep})" , 21, -4., 4.); lepWMass_ = fs->make<TH1F>("lepWMass", "M (W_{lep}) [GeV]" , 25, 0., 200.); lepTopPt_ = fs->make<TH1F>("lepTopPt" , "p_{T} (t_{lep}) [GeV]", 25, 0. , 500.); lepTopEta_ = fs->make<TH1F>("lepTopEta" , "#eta (t_{lep})" , 21, -4., 4.); lepTopMass_ = fs->make<TH1F>("lepTopMass", "M (t_{lep}) [GeV]" , 40, 0. , 400.); hadWPullPt_ = fs->make<TH1F>("hadWPullPt" , "(p_{T,rec}-p_{T,gen})/p_{T,gen} (W_{had})" , 40, -1., 1.); hadWPullEta_ = fs->make<TH1F>("hadWPullEta" , "(#eta_{rec}-#eta_{gen})/#eta_{gen} (W_{had})", 40, -1., 1.); hadWPullMass_ = fs->make<TH1F>("hadWPullMass", "(M_{rec}-M_{gen})/M_{gen} (W_{had})" , 40, -1., 1.); hadTopPullPt_ = fs->make<TH1F>("hadTopPullPt" , "(p_{T,rec}-p_{T,gen})/p_{T,gen} (t_{had})" , 40, -1., 1.); hadTopPullEta_ = fs->make<TH1F>("hadTopPullEta" , "(#eta_{rec}-#eta_{gen})/#eta_{gen} (t_{had})", 40, -1., 1.); hadTopPullMass_ = fs->make<TH1F>("hadTopPullMass", "(M_{rec}-M_{gen})/M_{gen} (t_{had})" , 40, -1., 1.); lepWPullPt_ = fs->make<TH1F>("lepWPullPt" , "(p_{T,rec}-p_{T,gen})/p_{T,gen} (W_{lep})" , 40, -1., 1.); lepWPullEta_ = fs->make<TH1F>("lepWPullEta" , "(#eta_{rec}-#eta_{gen})/#eta_{gen} (W_{lep})", 40, -1., 1.); lepWPullMass_ = fs->make<TH1F>("lepWPullMass", "(M_{rec}-M_{gen})/M_{gen} (W_{lep})" , 40, -1., 1.); lepTopPullPt_ = fs->make<TH1F>("lepTopPullPt" , "(p_{T,rec}-p_{T,gen})/p_{T,gen} (t_{lep})" , 40, -1., 1.); lepTopPullEta_ = fs->make<TH1F>("lepTopPullEta" , "(#eta_{rec}-#eta_{gen})/#eta_{gen} (t_{lep})", 40, -1., 1.); lepTopPullMass_ = fs->make<TH1F>("lepTopPullMass", "(M_{rec}-M_{gen})/M_{gen} (t_{lep})" , 40, -1., 1.); topPairMass_ = fs->make<TH1F>("topPairMass", "M (t#bar{t})", 36, 340., 940.); topPairPullMass_ = fs->make<TH1F>("topPairPullMass", "(M_{rec}-M_{gen})/M_{gen} (t#bar{t})", 40, -1., 1.); genMatchDr_ = fs->make<TH1F>("genMatchDr", "GenMatch #Sigma#DeltaR", 40, 0., 4.); kinFitProb_ = fs->make<TH1F>("kinFitProb", "KinFit probability" , 50, 0., 1.); genMatchDrVsHadTopPullMass_ = fs->make<TH2F>("genMatchDrVsHadTopPullMass", "GenMatch #Sigma #Delta R vs. (M_{rec}-M_{gen})/M_{gen} (t_{had}))", 40, -1., 1., 40, 0., 4.); kinFitProbVsHadTopPullMass_ = fs->make<TH2F>("kinFitProbVsHadTopPullMass", "KinFit probability vs. (M_{rec}-M_{gen})/M_{gen} (t_{had}))", 40, -1., 1., 20, 0., 1.); }
void HypothesisAnalyzer::endJob | ( | void | ) | [private, virtual] |
TH1F* HypothesisAnalyzer::genMatchDr_ [private] |
Definition at line 63 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH2F* HypothesisAnalyzer::genMatchDrVsHadTopPullMass_ [private] |
Definition at line 66 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadTopEta_ [private] |
Definition at line 37 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadTopMass_ [private] |
Definition at line 38 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadTopPt_ [private] |
Definition at line 36 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadTopPullEta_ [private] |
Definition at line 41 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadTopPullMass_ [private] |
Definition at line 42 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadTopPullPt_ [private] |
Definition at line 40 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadWEta_ [private] |
Definition at line 29 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadWMass_ [private] |
Definition at line 30 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadWPt_ [private] |
Definition at line 28 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadWPullEta_ [private] |
Definition at line 33 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadWPullMass_ [private] |
Definition at line 34 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::hadWPullPt_ [private] |
Definition at line 32 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
const std::string HypothesisAnalyzer::hypoClassKey_ [private] |
Definition at line 23 of file HypothesisAnalyzer.h.
Referenced by analyze().
TH1F* HypothesisAnalyzer::kinFitProb_ [private] |
Definition at line 64 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH2F* HypothesisAnalyzer::kinFitProbVsHadTopPullMass_ [private] |
Definition at line 67 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepTopEta_ [private] |
Definition at line 53 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepTopMass_ [private] |
Definition at line 54 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepTopPt_ [private] |
Definition at line 52 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepTopPullEta_ [private] |
Definition at line 60 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepTopPullMass_ [private] |
Definition at line 61 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepTopPullPt_ [private] |
Definition at line 59 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepWEta_ [private] |
Definition at line 45 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepWMass_ [private] |
Definition at line 46 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepWPt_ [private] |
Definition at line 44 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepWPullEta_ [private] |
Definition at line 49 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepWPullMass_ [private] |
Definition at line 50 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::lepWPullPt_ [private] |
Definition at line 48 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::neutrinoEta_ [private] |
Definition at line 25 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::neutrinoPullEta_ [private] |
Definition at line 26 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
const edm::InputTag HypothesisAnalyzer::semiLepEvt_ [private] |
Definition at line 22 of file HypothesisAnalyzer.h.
Referenced by analyze().
TH1F* HypothesisAnalyzer::topPairMass_ [private] |
Definition at line 56 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* HypothesisAnalyzer::topPairPullMass_ [private] |
Definition at line 57 of file HypothesisAnalyzer.h.
Referenced by analyze(), and beginJob().