ATLAS Offline Software
Loading...
Searching...
No Matches
FlavourUncertaintyComponent.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef JETUNCERTAINTIES_FLAVOURUNCERTAINTYCOMPONENT_H
6#define JETUNCERTAINTIES_FLAVOURUNCERTAINTYCOMPONENT_H
7
9
10namespace jet
11{
12
14{
15 public:
16 // Constructor/destructor
18 const TString& jetType,
19 const TString& analysisRootFileName,
20 const TString& defaultAnalysisRootFileName,
21 const TString& path,
22 const TString& calibArea,
23 const bool absEtaGluonFraction,
24 const TString& analysisHistPattern="",
25 const TString& NjetAccessorName="Njet"
26 );
28 virtual FlavourUncertaintyComponent* clone() const;
30 virtual StatusCode initialize(TFile* histFile);
31
32 // Extra information retrieval methods
34
35 protected:
36
37 // Uncertainty/validity retrieval helper methods
38 virtual bool getValidityImpl(const xAOD::Jet& jet, const xAOD::EventInfo& eInfo) const;
39 virtual double getUncertaintyImpl(const xAOD::Jet& jet, const xAOD::EventInfo& eInfo) const;
40
41 private:
42 // Default constructor is forbidden
43 FlavourUncertaintyComponent(const std::string& name = "");
44
46
47 // Additional private members
49 const TString m_jetType;
50 const TString m_analysisFileName;
51 const TString m_analysisHistPattern;
52 const TString m_defAnaFileName;
53 const TString m_path;
54 const TString m_calibArea;
55 const bool m_absEta;
57 const TString m_secondUncName;
58
60
61 // Large-R flags to only apply to specific jet types (as other jet types handled by topology uncertainties)
62 std::vector<LargeRJetTruthLabel::TypeEnum> m_largeRJetTruthLabels;
63
70
71 // Analysis histograms from analysis root file
72 std::vector<UncertaintyHistogram*> m_gluonFractionHists;
73 std::vector<UncertaintyHistogram*> m_gluonFractionErrorHists;
74
75 // Wrappers for special flavour histograms
76 double getFlavourResponseUncertainty(const xAOD::Jet& jet, const xAOD::EventInfo& eInfo) const;
77 double getFlavourCompositionUncertainty(const xAOD::Jet& jet, const xAOD::EventInfo& eInfo) const;
78 double getBJESUncertainty(const xAOD::Jet& jet, const xAOD::EventInfo& eInfo) const;
79
80 double getGluonFraction(const double pT, const double eta, const int nJets) const;
81 double getGluonFractionError(const double pT, const double eta, const int nJets) const;
82 double getGluonResponseDifference(const double pT, const double eta) const;
83 double getGluonResponseBaseline(const double pT, const double eta) const;
84 double getQuarkResponseBaseline(const double pT, const double eta) const;
85
86 // Private helper indices and functions
87 StatusCode readNjetsHistograms(std::vector<UncertaintyHistogram*>& hists, const std::vector<TString>& histKeys);
88 StatusCode getNjetFromKey(const TString& key, int& nJets) const;
89 StatusCode checkNjetsInput(int& nJets) const;
90 bool isBjet(const xAOD::Jet& jet) const;
91 void getGluonKeys(TFile* analysisFile, std::vector<TString>& gluonFractionKeys, std::vector<TString>& gluonFractionErrorKeys) const;
92};
93
94} // end jet namespace
95
96#endif
Scalar eta() const
pseudorapidity method
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
SG::AuxElement::Accessor< int > m_largeRJetTruthLabelAccessor
SG::AuxElement::Accessor< char > m_BjetAccessor
SG::AuxElement::Accessor< int > m_NjetAccessor
std::vector< UncertaintyHistogram * > m_gluonFractionErrorHists
double getGluonFractionError(const double pT, const double eta, const int nJets) const
StatusCode getNjetFromKey(const TString &key, int &nJets) const
double getBJESUncertainty(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
virtual FlavourUncertaintyComponent * clone() const
double getGluonFraction(const double pT, const double eta, const int nJets) const
double getGluonResponseBaseline(const double pT, const double eta) const
std::vector< UncertaintyHistogram * > m_gluonFractionHists
std::vector< LargeRJetTruthLabel::TypeEnum > m_largeRJetTruthLabels
void getGluonKeys(TFile *analysisFile, std::vector< TString > &gluonFractionKeys, std::vector< TString > &gluonFractionErrorKeys) const
virtual FlavourComp::TypeEnum getFlavourType() const
double getQuarkResponseBaseline(const double pT, const double eta) const
FlavourUncertaintyComponent(const ComponentHelper &component, const TString &jetType, const TString &analysisRootFileName, const TString &defaultAnalysisRootFileName, const TString &path, const TString &calibArea, const bool absEtaGluonFraction, const TString &analysisHistPattern="", const TString &NjetAccessorName="Njet")
StatusCode readNjetsHistograms(std::vector< UncertaintyHistogram * > &hists, const std::vector< TString > &histKeys)
double getFlavourResponseUncertainty(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
double getGluonResponseDifference(const double pT, const double eta) const
virtual double getUncertaintyImpl(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
virtual bool getValidityImpl(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
bool isBjet(const xAOD::Jet &jet) const
double getFlavourCompositionUncertainty(const xAOD::Jet &jet, const xAOD::EventInfo &eInfo) const
UncertaintyComponent(const ComponentHelper &component, const size_t numHist=1)
Jet_v1 Jet
Definition of the current "jet version".
EventInfo_v1 EventInfo
Definition of the latest event info version.
void initialize()