CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
HypothesisAnalyzer Class Reference

#include <HypothesisAnalyzer.h>

Inheritance diagram for HypothesisAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 HypothesisAnalyzer (const edm::ParameterSet &)
 
 ~HypothesisAnalyzer () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 

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::EDGetTokenT< TtSemiLeptonicEventsemiLepEvtToken_
 
TH1F * topPairMass_
 
TH1F * topPairPullMass_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 11 of file HypothesisAnalyzer.h.

Constructor & Destructor Documentation

HypothesisAnalyzer::HypothesisAnalyzer ( const edm::ParameterSet cfg)
explicit

Definition at line 9 of file HypothesisAnalyzer.cc.

10  : semiLepEvtToken_(consumes<TtSemiLeptonicEvent>(cfg.getParameter<edm::InputTag>("semiLepEvent"))),
11  hypoClassKey_(cfg.getParameter<std::string>("hypoClassKey")) {}
T getParameter(std::string const &) const
const std::string hypoClassKey_
const edm::EDGetTokenT< TtSemiLeptonicEvent > semiLepEvtToken_
HypothesisAnalyzer::~HypothesisAnalyzer ( )
inlineoverride

Definition at line 14 of file HypothesisAnalyzer.h.

References analyze(), beginJob(), and endJob().

14 {};

Member Function Documentation

void HypothesisAnalyzer::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
overrideprivate

Definition at line 13 of file HypothesisAnalyzer.cc.

References reco::Candidate::eta(), TtEvent::fitProb(), genMatchDr_, genMatchDrVsHadTopPullMass_, TtEvent::genMatchSumDR(), TtSemiLeptonicEvent::hadronicDecayTop(), TtSemiLeptonicEvent::hadronicDecayW(), hadTopEta_, hadTopMass_, hadTopPt_, hadTopPullEta_, hadTopPullMass_, hadTopPullPt_, hadWEta_, hadWMass_, hadWPt_, hadWPullEta_, hadWPullMass_, hadWPullPt_, hypoClassKey_, TtEvent::isHypoValid(), kinFitProb_, kinFitProbVsHadTopPullMass_, TtSemiLeptonicEvent::leptonicDecayTop(), TtSemiLeptonicEvent::leptonicDecayW(), lepTopEta_, lepTopMass_, lepTopPt_, lepTopPullEta_, lepTopPullMass_, lepTopPullPt_, lepWEta_, lepWMass_, lepWPt_, lepWPullEta_, lepWPullMass_, lepWPullPt_, reco::Candidate::mass(), neutrinoEta_, neutrinoPullEta_, reco::Candidate::pt(), semiLepEvtToken_, TtSemiLeptonicEvent::singleNeutrino(), TtEvent::topPair(), topPairMass_, and topPairPullMass_.

Referenced by ~HypothesisAnalyzer().

13  {
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 }
const reco::Candidate * leptonicDecayW(const std::string &key, const unsigned &cmb=0) const
get leptonic W of the given hypothesis
const std::string hypoClassKey_
const reco::Candidate * hadronicDecayW(const std::string &key, const unsigned &cmb=0) const
get hadronic W of the given hypothesis
const reco::Candidate * hadronicDecayTop(const std::string &key, const unsigned &cmb=0) const
get hadronic top of the given hypothesis
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double genMatchSumDR(const unsigned &cmb=0) const
return the sum dr of the generator match if available; -1 else
Definition: TtEvent.h:120
const edm::EDGetTokenT< TtSemiLeptonicEvent > semiLepEvtToken_
const reco::Candidate * leptonicDecayTop(const std::string &key, const unsigned &cmb=0) const
get leptonic top of the given hypothesis
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
virtual double mass() const =0
mass
double fitProb(const unsigned &cmb=0) const
return the fit probability of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:132
bool isHypoValid(const std::string &key, const unsigned &cmb=0) const
check if hypothesis &#39;cmb&#39; within the hypothesis class was valid; if not it lead to an empty Composite...
Definition: TtEvent.h:82
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:147
const reco::Candidate * singleNeutrino(const std::string &key, const unsigned &cmb=0) const
get leptonic light quark of the given hypothesis
void HypothesisAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 126 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_, TFileService::make(), neutrinoEta_, neutrinoPullEta_, topPairMass_, and topPairPullMass_.

Referenced by ~HypothesisAnalyzer().

126  {
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 }
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void HypothesisAnalyzer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 194 of file HypothesisAnalyzer.cc.

Referenced by ~HypothesisAnalyzer().

194 {}

Member Data Documentation

TH1F* HypothesisAnalyzer::genMatchDr_
private

Definition at line 62 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH2F* HypothesisAnalyzer::genMatchDrVsHadTopPullMass_
private

Definition at line 65 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadTopEta_
private

Definition at line 36 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadTopMass_
private

Definition at line 37 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadTopPt_
private

Definition at line 35 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadTopPullEta_
private

Definition at line 40 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadTopPullMass_
private

Definition at line 41 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadTopPullPt_
private

Definition at line 39 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadWEta_
private

Definition at line 28 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadWMass_
private

Definition at line 29 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadWPt_
private

Definition at line 27 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadWPullEta_
private

Definition at line 32 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadWPullMass_
private

Definition at line 33 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::hadWPullPt_
private

Definition at line 31 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

const std::string HypothesisAnalyzer::hypoClassKey_
private

Definition at line 22 of file HypothesisAnalyzer.h.

Referenced by analyze().

TH1F* HypothesisAnalyzer::kinFitProb_
private

Definition at line 63 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH2F* HypothesisAnalyzer::kinFitProbVsHadTopPullMass_
private

Definition at line 66 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepTopEta_
private

Definition at line 52 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepTopMass_
private

Definition at line 53 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepTopPt_
private

Definition at line 51 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepTopPullEta_
private

Definition at line 59 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepTopPullMass_
private

Definition at line 60 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepTopPullPt_
private

Definition at line 58 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepWEta_
private

Definition at line 44 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepWMass_
private

Definition at line 45 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepWPt_
private

Definition at line 43 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepWPullEta_
private

Definition at line 48 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepWPullMass_
private

Definition at line 49 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::lepWPullPt_
private

Definition at line 47 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::neutrinoEta_
private

Definition at line 24 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::neutrinoPullEta_
private

Definition at line 25 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

const edm::EDGetTokenT<TtSemiLeptonicEvent> HypothesisAnalyzer::semiLepEvtToken_
private

Definition at line 21 of file HypothesisAnalyzer.h.

Referenced by analyze().

TH1F* HypothesisAnalyzer::topPairMass_
private

Definition at line 55 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* HypothesisAnalyzer::topPairPullMass_
private

Definition at line 56 of file HypothesisAnalyzer.h.

Referenced by analyze(), and beginJob().