CMS 3D CMS Logo

HypothesisAnalyzer.cc
Go to the documentation of this file.
6 
8 
10  : semiLepEvtToken_(consumes<TtSemiLeptonicEvent>(cfg.getParameter<edm::InputTag>("semiLepEvent"))),
11  hypoClassKey_(cfg.getParameter<std::string>("hypoClassKey")) {}
12 
15  // get a handle for the TtSemiLeptonicEvent and a key to the hypothesis
17 
19  event.getByToken(semiLepEvtToken_, semiLepEvt);
20 
22  // check if hypothesis is available and valid in this event
24 
25  if (!semiLepEvt->isHypoValid(hypoClassKey_)) {
26  edm::LogInfo("HypothesisAnalyzer") << "Hypothesis " << hypoClassKey_ << " not valid for this event";
27  return;
28  }
29 
31  // get reconstructed top quarks, W bosons, the top pair and the neutrino from the hypothesis
33 
34  const reco::Candidate* topPair = semiLepEvt->topPair(hypoClassKey_);
35  const reco::Candidate* lepTop = semiLepEvt->leptonicDecayTop(hypoClassKey_);
36  const reco::Candidate* lepW = semiLepEvt->leptonicDecayW(hypoClassKey_);
37  const reco::Candidate* hadTop = semiLepEvt->hadronicDecayTop(hypoClassKey_);
38  const reco::Candidate* hadW = semiLepEvt->hadronicDecayW(hypoClassKey_);
39  const reco::Candidate* neutrino = semiLepEvt->singleNeutrino(hypoClassKey_);
40 
42  // fill simple histograms with kinematic variables of the reconstructed particles
44 
45  if (topPair)
46  topPairMass_->Fill(topPair->mass());
47  if (hadW) {
48  hadWPt_->Fill(hadW->pt());
49  hadWEta_->Fill(hadW->eta());
50  hadWMass_->Fill(hadW->mass());
51  }
52  if (hadTop) {
53  hadTopPt_->Fill(hadTop->pt());
54  hadTopEta_->Fill(hadTop->eta());
55  hadTopMass_->Fill(hadTop->mass());
56  }
57  if (lepW) {
58  lepWPt_->Fill(lepW->pt());
59  lepWEta_->Fill(lepW->eta());
60  lepWMass_->Fill(lepW->mass());
61  }
62  if (lepTop) {
63  lepTopPt_->Fill(lepTop->pt());
64  lepTopEta_->Fill(lepTop->eta());
65  lepTopMass_->Fill(lepTop->mass());
66  }
67  if (neutrino)
68  neutrinoEta_->Fill(neutrino->eta());
69 
71  // get corresponding genParticles
73 
74  const math::XYZTLorentzVector* genTopPair = semiLepEvt->topPair();
75  const reco::Candidate* genHadTop = semiLepEvt->hadronicDecayTop();
76  const reco::Candidate* genHadW = semiLepEvt->hadronicDecayW();
77  const reco::Candidate* genLepTop = semiLepEvt->leptonicDecayTop();
78  const reco::Candidate* genLepW = semiLepEvt->leptonicDecayW();
79  const reco::Candidate* genNeutrino = semiLepEvt->singleNeutrino();
80 
82  // fill pull histograms of kinematic variables with respect to the generated particles
84 
85  if (topPair && genTopPair)
86  topPairPullMass_->Fill((topPair->mass() - genTopPair->mass()) / genTopPair->mass());
87  if (hadW && genHadW) {
88  hadWPullPt_->Fill((hadW->pt() - genHadW->pt()) / genHadW->pt());
89  hadWPullEta_->Fill((hadW->eta() - genHadW->eta()) / genHadW->eta());
90  hadWPullMass_->Fill((hadW->mass() - genHadW->mass()) / genHadW->mass());
91  }
92 
93  if (hadTop && genHadTop) {
94  hadTopPullPt_->Fill((hadTop->pt() - genHadTop->pt()) / genHadTop->pt());
95  hadTopPullEta_->Fill((hadTop->eta() - genHadTop->eta()) / genHadTop->eta());
96  hadTopPullMass_->Fill((hadTop->mass() - genHadTop->mass()) / genHadTop->mass());
97  }
98  if (lepW && genLepW) {
99  lepWPullPt_->Fill((lepW->pt() - genLepW->pt()) / genLepW->pt());
100  lepWPullEta_->Fill((lepW->eta() - genLepW->eta()) / genLepW->eta());
101  lepWPullMass_->Fill((lepW->mass() - genLepW->mass()) / genLepW->mass());
102  }
103 
104  if (lepTop && genLepTop) {
105  lepTopPullPt_->Fill((lepTop->pt() - genLepTop->pt()) / genLepTop->pt());
106  lepTopPullEta_->Fill((lepTop->eta() - genLepTop->eta()) / genLepTop->eta());
107  lepTopPullMass_->Fill((lepTop->mass() - genLepTop->mass()) / genLepTop->mass());
108  }
109  if (neutrino && genNeutrino)
110  neutrinoPullEta_->Fill((neutrino->eta() - genNeutrino->eta()) / genNeutrino->eta());
111 
113  // fill histograms with variables describing the quality of the hypotheses
115 
116  genMatchDr_->Fill(semiLepEvt->genMatchSumDR());
117  kinFitProb_->Fill(semiLepEvt->fitProb());
118 
119  if (hadTop && genHadTop) {
120  genMatchDrVsHadTopPullMass_->Fill((hadTop->mass() - genHadTop->mass()) / genHadTop->mass(),
121  semiLepEvt->genMatchSumDR());
122  kinFitProbVsHadTopPullMass_->Fill((hadTop->mass() - genHadTop->mass()) / genHadTop->mass(), semiLepEvt->fitProb());
123  }
124 }
125 
128  if (!fs)
129  throw edm::Exception(edm::errors::Configuration, "TFile Service is not registered in cfg file");
130 
132  // book histograms
134 
135  neutrinoEta_ = fs->make<TH1F>("neutrinoEta", "#eta (neutrino)", 21, -4., 4.);
136  neutrinoPullEta_ = fs->make<TH1F>("neutrinoPullEta", "(#eta_{rec}-#eta_{gen})/#eta_{gen} (neutrino)", 40, -1., 1.);
137 
138  hadWPt_ = fs->make<TH1F>("hadWPt", "p_{T} (W_{had}) [GeV]", 25, 0., 500.);
139  hadWEta_ = fs->make<TH1F>("hadWEta", "#eta (W_{had})", 21, -4., 4.);
140  hadWMass_ = fs->make<TH1F>("hadWMass", "M (W_{had}) [GeV]", 25, 0., 200.);
141 
142  hadTopPt_ = fs->make<TH1F>("hadTopPt", "p_{T} (t_{had}) [GeV]", 25, 0., 500.);
143  hadTopEta_ = fs->make<TH1F>("hadTopEta", "#eta (t_{had})", 21, -4., 4.);
144  hadTopMass_ = fs->make<TH1F>("hadTopMass", "M (t_{had}) [GeV]", 40, 0., 400.);
145 
146  lepWPt_ = fs->make<TH1F>("lepWPt", "p_{t} (W_{lep}) [GeV]", 25, 0., 500.);
147  lepWEta_ = fs->make<TH1F>("lepWEta", "#eta (W_{lep})", 21, -4., 4.);
148  lepWMass_ = fs->make<TH1F>("lepWMass", "M (W_{lep}) [GeV]", 25, 0., 200.);
149 
150  lepTopPt_ = fs->make<TH1F>("lepTopPt", "p_{T} (t_{lep}) [GeV]", 25, 0., 500.);
151  lepTopEta_ = fs->make<TH1F>("lepTopEta", "#eta (t_{lep})", 21, -4., 4.);
152  lepTopMass_ = fs->make<TH1F>("lepTopMass", "M (t_{lep}) [GeV]", 40, 0., 400.);
153 
154  hadWPullPt_ = fs->make<TH1F>("hadWPullPt", "(p_{T,rec}-p_{T,gen})/p_{T,gen} (W_{had})", 40, -1., 1.);
155  hadWPullEta_ = fs->make<TH1F>("hadWPullEta", "(#eta_{rec}-#eta_{gen})/#eta_{gen} (W_{had})", 40, -1., 1.);
156  hadWPullMass_ = fs->make<TH1F>("hadWPullMass", "(M_{rec}-M_{gen})/M_{gen} (W_{had})", 40, -1., 1.);
157 
158  hadTopPullPt_ = fs->make<TH1F>("hadTopPullPt", "(p_{T,rec}-p_{T,gen})/p_{T,gen} (t_{had})", 40, -1., 1.);
159  hadTopPullEta_ = fs->make<TH1F>("hadTopPullEta", "(#eta_{rec}-#eta_{gen})/#eta_{gen} (t_{had})", 40, -1., 1.);
160  hadTopPullMass_ = fs->make<TH1F>("hadTopPullMass", "(M_{rec}-M_{gen})/M_{gen} (t_{had})", 40, -1., 1.);
161 
162  lepWPullPt_ = fs->make<TH1F>("lepWPullPt", "(p_{T,rec}-p_{T,gen})/p_{T,gen} (W_{lep})", 40, -1., 1.);
163  lepWPullEta_ = fs->make<TH1F>("lepWPullEta", "(#eta_{rec}-#eta_{gen})/#eta_{gen} (W_{lep})", 40, -1., 1.);
164  lepWPullMass_ = fs->make<TH1F>("lepWPullMass", "(M_{rec}-M_{gen})/M_{gen} (W_{lep})", 40, -1., 1.);
165 
166  lepTopPullPt_ = fs->make<TH1F>("lepTopPullPt", "(p_{T,rec}-p_{T,gen})/p_{T,gen} (t_{lep})", 40, -1., 1.);
167  lepTopPullEta_ = fs->make<TH1F>("lepTopPullEta", "(#eta_{rec}-#eta_{gen})/#eta_{gen} (t_{lep})", 40, -1., 1.);
168  lepTopPullMass_ = fs->make<TH1F>("lepTopPullMass", "(M_{rec}-M_{gen})/M_{gen} (t_{lep})", 40, -1., 1.);
169 
170  topPairMass_ = fs->make<TH1F>("topPairMass", "M (t#bar{t})", 36, 340., 940.);
171  topPairPullMass_ = fs->make<TH1F>("topPairPullMass", "(M_{rec}-M_{gen})/M_{gen} (t#bar{t})", 40, -1., 1.);
172 
173  genMatchDr_ = fs->make<TH1F>("genMatchDr", "GenMatch #Sigma#DeltaR", 40, 0., 4.);
174  kinFitProb_ = fs->make<TH1F>("kinFitProb", "KinFit probability", 50, 0., 1.);
175 
176  genMatchDrVsHadTopPullMass_ = fs->make<TH2F>("genMatchDrVsHadTopPullMass",
177  "GenMatch #Sigma #Delta R vs. (M_{rec}-M_{gen})/M_{gen} (t_{had}))",
178  40,
179  -1.,
180  1.,
181  40,
182  0.,
183  4.);
184  kinFitProbVsHadTopPullMass_ = fs->make<TH2F>("kinFitProbVsHadTopPullMass",
185  "KinFit probability vs. (M_{rec}-M_{gen})/M_{gen} (t_{had}))",
186  40,
187  -1.,
188  1.,
189  20,
190  0.,
191  1.);
192 }
193 
HypothesisAnalyzer::lepTopPullMass_
TH1F * lepTopPullMass_
Definition: HypothesisAnalyzer.h:60
HypothesisAnalyzer::hadWPullPt_
TH1F * hadWPullPt_
Definition: HypothesisAnalyzer.h:31
HypothesisAnalyzer::genMatchDrVsHadTopPullMass_
TH2F * genMatchDrVsHadTopPullMass_
Definition: HypothesisAnalyzer.h:65
HypothesisAnalyzer::hadTopPullMass_
TH1F * hadTopPullMass_
Definition: HypothesisAnalyzer.h:41
MessageLogger.h
HypothesisAnalyzer::hadTopPullPt_
TH1F * hadTopPullPt_
Definition: HypothesisAnalyzer.h:39
TtSemiLeptonicEvent::leptonicDecayTop
const reco::Candidate * leptonicDecayTop(const std::string &key, const unsigned &cmb=0) const
get leptonic top of the given hypothesis
Definition: TtSemiLeptonicEvent.h:72
TtEvent::isHypoValid
bool isHypoValid(const std::string &key, const unsigned &cmb=0) const
check if hypothesis 'cmb' within the hypothesis class was valid; if not it lead to an empty Composite...
Definition: TtEvent.h:78
HypothesisAnalyzer::semiLepEvtToken_
const edm::EDGetTokenT< TtSemiLeptonicEvent > semiLepEvtToken_
Definition: HypothesisAnalyzer.h:21
HypothesisAnalyzer::neutrinoEta_
TH1F * neutrinoEta_
Definition: HypothesisAnalyzer.h:24
reco::Candidate::mass
virtual double mass() const =0
mass
reco::Candidate::eta
virtual double eta() const =0
momentum pseudorapidity
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
HypothesisAnalyzer::lepTopPt_
TH1F * lepTopPt_
Definition: HypothesisAnalyzer.h:51
reco::Candidate::pt
virtual double pt() const =0
transverse momentum
TtEvent::genMatchSumDR
double genMatchSumDR(const unsigned &cmb=0) const
return the sum dr of the generator match if available; -1 else
Definition: TtEvent.h:116
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
TtEvent::fitProb
double fitProb(const unsigned &cmb=0) const
return the fit probability of hypothesis 'cmb' if available; -1 else
Definition: TtEvent.h:128
edm::Handle
Definition: AssociativeIterator.h:50
singleTopDQM_cfi.setup
setup
Definition: singleTopDQM_cfi.py:37
HypothesisAnalyzer::endJob
void endJob() override
Definition: HypothesisAnalyzer.cc:194
HypothesisAnalyzer::hadWPt_
TH1F * hadWPt_
Definition: HypothesisAnalyzer.h:27
TtSemiLeptonicEvent
Class derived from the TtEvent for the semileptonic decay channel.
Definition: TtSemiLeptonicEvent.h:24
HypothesisAnalyzer::lepTopPullEta_
TH1F * lepTopPullEta_
Definition: HypothesisAnalyzer.h:59
HypothesisAnalyzer::lepTopPullPt_
TH1F * lepTopPullPt_
Definition: HypothesisAnalyzer.h:58
HypothesisAnalyzer::lepTopMass_
TH1F * lepTopMass_
Definition: HypothesisAnalyzer.h:53
Service.h
HypothesisAnalyzer::genMatchDr_
TH1F * genMatchDr_
Definition: HypothesisAnalyzer.h:62
TtSemiLeptonicEvent::hadronicDecayTop
const reco::Candidate * hadronicDecayTop(const std::string &key, const unsigned &cmb=0) const
get hadronic top of the given hypothesis
Definition: TtSemiLeptonicEvent.h:32
HypothesisAnalyzer.h
HypothesisAnalyzer::beginJob
void beginJob() override
Definition: HypothesisAnalyzer.cc:126
TtSemiLeptonicEvent::hadronicDecayW
const reco::Candidate * hadronicDecayW(const std::string &key, const unsigned &cmb=0) const
get hadronic W of the given hypothesis
Definition: TtSemiLeptonicEvent.h:48
HypothesisAnalyzer::hadWPullMass_
TH1F * hadWPullMass_
Definition: HypothesisAnalyzer.h:33
HypothesisAnalyzer::lepWPullPt_
TH1F * lepWPullPt_
Definition: HypothesisAnalyzer.h:47
HypothesisAnalyzer::lepTopEta_
TH1F * lepTopEta_
Definition: HypothesisAnalyzer.h:52
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
TFileService.h
HypothesisAnalyzer::hadTopMass_
TH1F * hadTopMass_
Definition: HypothesisAnalyzer.h:37
edm::ParameterSet
Definition: ParameterSet.h:47
HypothesisAnalyzer::hadWEta_
TH1F * hadWEta_
Definition: HypothesisAnalyzer.h:28
Event.h
HypothesisAnalyzer::hadTopPt_
TH1F * hadTopPt_
Definition: HypothesisAnalyzer.h:35
HypothesisAnalyzer::topPairMass_
TH1F * topPairMass_
Definition: HypothesisAnalyzer.h:55
HypothesisAnalyzer::lepWMass_
TH1F * lepWMass_
Definition: HypothesisAnalyzer.h:45
HypothesisAnalyzer::neutrinoPullEta_
TH1F * neutrinoPullEta_
Definition: HypothesisAnalyzer.h:25
HypothesisAnalyzer::hadTopPullEta_
TH1F * hadTopPullEta_
Definition: HypothesisAnalyzer.h:40
edm::Service< TFileService >
edm::EventSetup
Definition: EventSetup.h:57
HypothesisAnalyzer::lepWPullMass_
TH1F * lepWPullMass_
Definition: HypothesisAnalyzer.h:49
TtSemiLeptonicEvent::leptonicDecayW
const reco::Candidate * leptonicDecayW(const std::string &key, const unsigned &cmb=0) const
get leptonic W of the given hypothesis
Definition: TtSemiLeptonicEvent.h:88
HypothesisAnalyzer::kinFitProbVsHadTopPullMass_
TH2F * kinFitProbVsHadTopPullMass_
Definition: HypothesisAnalyzer.h:66
TtSemiLeptonicEvent::singleNeutrino
const reco::Candidate * singleNeutrino(const std::string &key, const unsigned &cmb=0) const
get leptonic light quark of the given hypothesis
Definition: TtSemiLeptonicEvent.h:96
looper.cfg
cfg
Definition: looper.py:297
reco::Candidate
Definition: Candidate.h:27
HypothesisAnalyzer::topPairPullMass_
TH1F * topPairPullMass_
Definition: HypothesisAnalyzer.h:56
std
Definition: JetResolutionObject.h:76
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Exception
Definition: hltDiff.cc:246
HypothesisAnalyzer::lepWEta_
TH1F * lepWEta_
Definition: HypothesisAnalyzer.h:44
HypothesisAnalyzer::hadTopEta_
TH1F * hadTopEta_
Definition: HypothesisAnalyzer.h:36
HypothesisAnalyzer::HypothesisAnalyzer
HypothesisAnalyzer(const edm::ParameterSet &)
Definition: HypothesisAnalyzer.cc:9
HypothesisAnalyzer::hadWMass_
TH1F * hadWMass_
Definition: HypothesisAnalyzer.h:29
HypothesisAnalyzer::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: HypothesisAnalyzer.cc:13
ParameterSet.h
event
Definition: event.py:1
HypothesisAnalyzer::lepWPullEta_
TH1F * lepWPullEta_
Definition: HypothesisAnalyzer.h:48
edm::Event
Definition: Event.h:73
edm::errors::Configuration
Definition: EDMException.h:36
TtEvent::topPair
const reco::Candidate * topPair(const std::string &key, const unsigned &cmb=0) const
get combined 4-vector of top and topBar of the given hypothesis
Definition: TtEvent.h:143
HypothesisAnalyzer::lepWPt_
TH1F * lepWPt_
Definition: HypothesisAnalyzer.h:43
HypothesisAnalyzer::hypoClassKey_
const std::string hypoClassKey_
Definition: HypothesisAnalyzer.h:22
TFileService::make
T * make(const Args &... args) const
make new ROOT object
Definition: TFileService.h:64
HypothesisAnalyzer::kinFitProb_
TH1F * kinFitProb_
Definition: HypothesisAnalyzer.h:63
HypothesisAnalyzer::hadWPullEta_
TH1F * hadWPullEta_
Definition: HypothesisAnalyzer.h:32