ATLAS Offline Software
Public Member Functions | Private Attributes | List of all members
top::JetResponsePlots Class Reference

An example of how to quickly make some plots at a certain point in the cutflow. More...

#include <JetResponsePlots.h>

Inheritance diagram for top::JetResponsePlots:
Collaboration diagram for top::JetResponsePlots:

Public Member Functions

 JetResponsePlots (const std::string &name, TFile *outputFile, const std::string &params, std::shared_ptr< top::TopConfig > config, EL::Worker *wk=nullptr)
 Setup some example plots. More...
 
bool apply (const top::Event &event) const override
 Fill the histograms. More...
 
void FillHistograms (const double w_event, const top::Event &event) const
 Helper function to fill the histograms. More...
 
std::string name () const override
 Return the name for the cutflow table. More...
 
virtual bool applyParticleLevel (const top::ParticleLevelEvent &) const
 This does stuff based on the information in a particle level event. More...
 

Private Attributes

std::shared_ptr< PlotManagerm_hists
 
std::size_t m_nominalHashValue
 
float m_deltaR
 
int m_bins
 
float m_min
 
float m_max
 
std::vector< double > m_ptBinning
 
std::shared_ptr< top::TopConfigm_config
 
ToolHandle< PMGTools::IPMGTruthWeightToolm_PMGTruthWeights
 

Detailed Description

An example of how to quickly make some plots at a certain point in the cutflow.

Definition at line 31 of file JetResponsePlots.h.

Constructor & Destructor Documentation

◆ JetResponsePlots()

top::JetResponsePlots::JetResponsePlots ( const std::string &  name,
TFile *  outputFile,
const std::string &  params,
std::shared_ptr< top::TopConfig config,
EL::Worker wk = nullptr 
)

Setup some example plots.

Add a bunch of histograms.

Parameters
nameThe name of the directory to store histograms in, in the output file. e.g. you might have ee, mumu and emu.
outputFileThe output file. Needs setting up at the very start so that we can attach the files.
paramsThe arguments, e.g. for the binning of the plots.
configInstance of TopConfig
wkOnly used by EventLoop, ok as nullptr as default.

Definition at line 24 of file JetResponsePlots.cxx.

28  :
30  m_deltaR(0.3),
31  m_bins(30),
32  m_min(-3),
33  m_max(3),
35  m_PMGTruthWeights(nullptr) {
36 
38  m_nominalHashValue = nominal.hash();
39 
40  //retrieve PMGTruthWeights
41  static const std::string truthWeightToolName = "PMGTruthWeightTool";
42  if (asg::ToolStore::contains<PMGTools::IPMGTruthWeightTool>(truthWeightToolName)) m_PMGTruthWeights =
43  asg::ToolStore::get<PMGTools::IPMGTruthWeightTool>(truthWeightToolName);
44  else {
45  std::unique_ptr<PMGTools::PMGTruthWeightTool> pmgtruthtool = std::make_unique<PMGTools::PMGTruthWeightTool>(truthWeightToolName);
46  top::check(pmgtruthtool->initialize(), "Failed to initialize " + truthWeightToolName);
47  m_PMGTruthWeights = pmgtruthtool.release();
48  }
49  //decoding arguments
50  std::istringstream stream(params);
51  std::string s;
52  while (stream >> s) {
53  if (s.substr(0,7) == "deltaR:") {
54  m_deltaR = std::stof(s.substr(7, s.size() - 7));
55  m_config->setJetResponseMatchingDeltaR(m_deltaR);
56  } else if (s.substr(0,5) == "bins:") {
57  m_bins = std::stoi(s.substr(5, s.size() - 5));
58  } else if (s.substr(0,4) == "min:") {
59  m_min = std::stof(s.substr(4, s.size() - 4));
60  } else if (s.substr(0,4) == "max:") {
61  m_max = std::stof(s.substr(4, s.size() - 4));
62  } else if (s.substr(0,10) == "ptBinning:") {
63  const std::string tmp = s.substr(10, s.size() - 4);
64  // split by the delimiter
65  std::vector<std::string> split;
66  tokenize(tmp, split, ",");
67  for (const std::string& istring : split) {
68  m_ptBinning.emplace_back(std::stof(istring));
69  }
70  } else {
71  throw std::runtime_error{"ERROR: Can't understand argument " + s + "For JETRESPONSEPLOTS"};
72  }
73  }
74 
75  if (m_ptBinning.size() < 2) {
76  throw std::invalid_argument{"pT binning (x axis) size is smaller than 2"};
77  }
78 
79  m_hists = std::make_shared<PlotManager>(name + "/JetResponsePlots", outputFile, wk);
80  m_hists->addHist("JetResponse", ";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
81  m_hists->addHist("JetResponse_bb", ";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
82  m_hists->addHist("JetResponse_b", ";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
83  m_hists->addHist("JetResponse_c", ";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
84  m_hists->addHist("JetResponse_tau", ";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
85  m_hists->addHist("JetResponse_q", ";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
86  m_hists->addHist("JetResponse_g", ";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
87  m_hists->addHist("JetResponse_other",";p_{T}^{truth} [GeV];p_{T}^{reco}/p_{T}^{truth};Events", m_ptBinning.size() - 1, m_ptBinning.data(), m_bins, m_min, m_max);
88  }

Member Function Documentation

◆ apply()

bool top::JetResponsePlots::apply ( const top::Event event) const
overridevirtual

Fill the histograms.

Returns
True because it doesn't select any events.

Implements top::EventSelectorBase.

Definition at line 90 of file JetResponsePlots.cxx.

90  {
91  //only MC
92  if (!top::isSimulation(event)) return true;
93 
94  // only nominal
95  if (event.m_hashValue != m_nominalHashValue) return true;
96 
97  // do loose or tight events only if requested
98  if (event.m_isLoose && !m_config->doLooseEvents()) return true;
99 
100  if (!event.m_isLoose && !m_config->doTightEvents()) return true;
101 
102  const double nominalWeight = event.m_info->auxdata<float>("AnalysisTop_eventWeight");
103 
104  FillHistograms(nominalWeight, event);
105 
106  return true;
107  }

◆ applyParticleLevel()

virtual bool top::EventSelectorBase::applyParticleLevel ( const top::ParticleLevelEvent ) const
inlinevirtualinherited

This does stuff based on the information in a particle level event.

The idea is that you implement this to return either true or false, based on the information held within the top::ParticleLevelEvent. If this function returns true, then the event is kept, otherwise it is removed. The function has a default implementation (which returns true) because it is expected that many EventSelector objects do not operate on ParticleLevelEvent objects.

Parameters
top::ParticleLevelEventthe current particle level event.
trueif the event should be kept (i.e. it passed the selector criteria), false otherwise.

Reimplemented in top::JetNGhostSelector, top::PrintEventSelector, top::PseudoTopRecoRun, top::NElectronNMuonTightSelector, top::NElectronNMuonSelector, top::NFwdElectronSelector, top::HTSelector, top::OSLeptonTightSelector, top::MLLSelector, top::MWTSelector, top::NElectronTightSelector, top::NFwdElectronTightSelector, top::NMuonTightSelector, top::OSLeptonSelector, top::METMWTSelector, top::METSelector, top::MLLWindow, top::NElectronSelector, top::NJetSelector, top::NMuonSelector, top::NPhotonSelector, top::NSoftMuonSelector, top::NTauSelector, top::SSLeptonTightSelector, top::SSLeptonSelector, top::ParticleLevelSelector, top::RecoLevelSelector, top::NVarRCJetSelector, top::NLargeJetSelector, and top::NRCJetSelector.

Definition at line 73 of file EventSelectorBase.h.

73 {return true;}

◆ FillHistograms()

void top::JetResponsePlots::FillHistograms ( const double  w_event,
const top::Event event 
) const

Helper function to fill the histograms.

Definition at line 109 of file JetResponsePlots.cxx.

110  {
111 
112  for (const auto* const jetPtr : event.m_jets) {
113  int jet_flavor = -99;
114  bool status = jetPtr->getAttribute<int>("HadronConeExclExtendedTruthLabelID", jet_flavor);
115  if (!status) continue;
116 
117  float matchedPt = -9999;
118  status = jetPtr->getAttribute<float>("AnalysisTop_MatchedTruthJetPt", matchedPt);
119  if (!status) continue;
120 
121  // protection against division by zero
122  if (matchedPt < 1e-6) continue;
123 
124  const float response = jetPtr->pt() / matchedPt;
125 
126  static_cast<TH2D*>(m_hists->hist("JetResponse"))->Fill(matchedPt/1e3, response, w_event);
127  if (jet_flavor == 55) {
128  static_cast<TH2D*>(m_hists->hist("JetResponse_bb"))->Fill(matchedPt/1e3, response, w_event);
129  } else if (jet_flavor == 5 || jet_flavor == 54) {
130  static_cast<TH2D*>(m_hists->hist("JetResponse_b"))->Fill(matchedPt/1e3, response, w_event);
131  } else if (jet_flavor == 4 || jet_flavor == 44) {
132  static_cast<TH2D*>(m_hists->hist("JetResponse_c"))->Fill(matchedPt/1e3, response, w_event);
133  } else if (jet_flavor == 15) {
134  static_cast<TH2D*>(m_hists->hist("JetResponse_tau"))->Fill(matchedPt/1e3, response, w_event);
135  } else {
136  status = jetPtr->getAttribute<int>("PartonTruthLabelID", jet_flavor);
137  if (!status) continue;
138  if (jet_flavor >= 1 && jet_flavor <= 3) {
139  static_cast<TH2D*>(m_hists->hist("JetResponse_q"))->Fill(matchedPt/1e3, response, w_event);
140  } else if (jet_flavor == 21) {
141  static_cast<TH2D*>(m_hists->hist("JetResponse_g"))->Fill(matchedPt/1e3, response, w_event);
142  } else {
143  static_cast<TH2D*>(m_hists->hist("JetResponse_other"))->Fill(matchedPt/1e3, response, w_event);
144  }
145  }
146  }
147  }

◆ name()

std::string top::JetResponsePlots::name ( ) const
overridevirtual

Return the name for the cutflow table.

Returns
The word JETRESPONSEPLOTS.

Implements top::EventSelectorBase.

Definition at line 149 of file JetResponsePlots.cxx.

149  {
150  return "JETRESPONSEPLOTS";
151  }

Member Data Documentation

◆ m_bins

int top::JetResponsePlots::m_bins
private

Definition at line 82 of file JetResponsePlots.h.

◆ m_config

std::shared_ptr<top::TopConfig> top::JetResponsePlots::m_config
private

Definition at line 88 of file JetResponsePlots.h.

◆ m_deltaR

float top::JetResponsePlots::m_deltaR
private

Definition at line 81 of file JetResponsePlots.h.

◆ m_hists

std::shared_ptr<PlotManager> top::JetResponsePlots::m_hists
private

Definition at line 75 of file JetResponsePlots.h.

◆ m_max

float top::JetResponsePlots::m_max
private

Definition at line 84 of file JetResponsePlots.h.

◆ m_min

float top::JetResponsePlots::m_min
private

Definition at line 83 of file JetResponsePlots.h.

◆ m_nominalHashValue

std::size_t top::JetResponsePlots::m_nominalHashValue
private

Definition at line 78 of file JetResponsePlots.h.

◆ m_PMGTruthWeights

ToolHandle<PMGTools::IPMGTruthWeightTool> top::JetResponsePlots::m_PMGTruthWeights
private

Definition at line 91 of file JetResponsePlots.h.

◆ m_ptBinning

std::vector<double> top::JetResponsePlots::m_ptBinning
private

Definition at line 85 of file JetResponsePlots.h.


The documentation for this class was generated from the following files:
top::JetResponsePlots::m_PMGTruthWeights
ToolHandle< PMGTools::IPMGTruthWeightTool > m_PMGTruthWeights
Definition: JetResponsePlots.h:91
python.SystemOfUnits.s
int s
Definition: SystemOfUnits.py:131
top::JetResponsePlots::m_deltaR
float m_deltaR
Definition: JetResponsePlots.h:81
top::JetResponsePlots::m_hists
std::shared_ptr< PlotManager > m_hists
Definition: JetResponsePlots.h:75
top::JetResponsePlots::FillHistograms
void FillHistograms(const double w_event, const top::Event &event) const
Helper function to fill the histograms.
Definition: JetResponsePlots.cxx:109
top::JetResponsePlots::m_bins
int m_bins
Definition: JetResponsePlots.h:82
response
MDT_Response response
Definition: MDT_ResponseTest.cxx:28
top::tokenize
void tokenize(const std::string &input, Container &output, const std::string &delimiters=" ", bool trim_empty=false)
Tokenize an input string using a set of delimiters.
Definition: Tokenize.h:24
CP::SystematicSet
Class to wrap a set of SystematicVariations.
Definition: SystematicSet.h:31
PMGTools::PMGTruthWeightTool::initialize
virtual StatusCode initialize() override
Function initialising the tool.
Definition: PMGTruthWeightTool.cxx:31
top::JetResponsePlots::m_config
std::shared_ptr< top::TopConfig > m_config
Definition: JetResponsePlots.h:88
AthenaPoolTestWrite.stream
string stream
Definition: AthenaPoolTestWrite.py:12
top::JetResponsePlots::name
std::string name() const override
Return the name for the cutflow table.
Definition: JetResponsePlots.cxx:149
config
Definition: PhysicsAnalysis/AnalysisCommon/AssociationUtils/python/config.py:1
compareGeometries.outputFile
string outputFile
Definition: compareGeometries.py:25
CheckAppliedSFs.e3
e3
Definition: CheckAppliedSFs.py:264
event
POOL::TEvent event(POOL::TEvent::kClassAccess)
top::nominal
@ nominal
Definition: ScaleFactorRetriever.h:29
top::check
void check(bool thingToCheck, const std::string &usefulFailureMessage)
Print an error message and terminate if thingToCheck is false.
Definition: EventTools.cxx:15
top::JetResponsePlots::m_ptBinning
std::vector< double > m_ptBinning
Definition: JetResponsePlots.h:85
TH2D
Definition: rootspy.cxx:430
top::isSimulation
bool isSimulation(const top::Event &event)
Is this event MC simulation (True) or data (False)?
Definition: EventTools.cxx:52
DeMoUpdate.tmp
string tmp
Definition: DeMoUpdate.py:1167
top::JetResponsePlots::m_min
float m_min
Definition: JetResponsePlots.h:83
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
top::JetResponsePlots::m_max
float m_max
Definition: JetResponsePlots.h:84
top::JetResponsePlots::m_nominalHashValue
std::size_t m_nominalHashValue
Definition: JetResponsePlots.h:78
merge.status
status
Definition: merge.py:17
PowhegControl_ttFCNC_NLO.params
params
Definition: PowhegControl_ttFCNC_NLO.py:226
Trk::split
@ split
Definition: LayerMaterialProperties.h:38