ATLAS Offline Software
Loading...
Searching...
No Matches
ResoTriggerMuonPlots Class Reference

#include <ResoTriggerMuonPlots.h>

Inheritance diagram for ResoTriggerMuonPlots:
Collaboration diagram for ResoTriggerMuonPlots:

Public Member Functions

 ResoTriggerMuonPlots (PlotBase *pParent, const std::string &sDir, std::string sType="")
void fill (const xAOD::Muon &Trigmu, const xAOD::Muon &Recomu)
void fill (const xAOD::L2StandAloneMuon &L2SAmu, const xAOD::Muon &Recomu)
void fill (const xAOD::L2CombinedMuon &L2CBmu, const xAOD::Muon &Recomu)
void initialize ()
void finalize ()
void setDetailLevel (int iDetailLevel)
void RegisterSubPlot (PlotBase *pPlotBase)
std::vector< HistDataretrieveBookedHistograms ()
 Retrieve all booked histograms.
std::vector< TreeDataretrieveBookedTrees ()
 Retrieve all booked trees.
std::vector< EfficiencyDataretrieveBookedEfficiencies ()
 Retrieve all booked efficiency objects.
TTree * BookTree (const std::string &name, bool prependDir=true)
 Book a TTree.
const std::string & getDirectory ()
Methods to book monitoring histograms

Note: methods starting with capitals should be deprecated in favour of camel-cased methods

TH1D * Book1D (const std::string &name, const std::string &labels, int nBins, float start, float end, bool prependDir=true)
 Book a TH1D histogram.
TH1D * Book1D (const std::string &name, TH1 *refHist, const std::string &labels, bool prependDir=true)
 Book a TH1D histogram using refHist as reference for number of bins and axis range.
TH2F * Book2D (const std::string &name, const std::string &labels, int nBinsX, float startX, float endX, int nBinsY, float startY, float endY, bool prependDir=true)
 Book a TH2F histogram.
TH2F * Book2D (const std::string &name, TH2 *refHist, const std::string &labels, bool prependDir=true)
 Book a TH2D histogram using refHist as reference for number of bins and axis range.
TH2F * Book2D (const std::string &name, const std::string &labels, int nBinsX, Double_t *binsX, int nBinsY, Double_t startY, Double_t endY, bool prependDir=true)
 Book a TH2F histogram with variable x axis binning.
TH3F * Book3D (const std::string &name, const std::string &labels, int nBinsX, float startX, float endX, int nBinsY, float startY, float endY, int nBinsZ, float startZ, float endZ, bool prependDir=true)
 Book a TH3F histogram.
TH3F * Book3D (const std::string &name, TH3 *refHist, const std::string &labels, bool prependDir=true)
 Book a TH3F histogram using refHist as reference for number of bins and axis range.
TProfile * BookTProfile (const std::string &name, const std::string &labels, int nBinsX, float startX, float endX, float startY=-1, float endY=-1, bool prependDir=true, bool useRMS=false)
 Book a TProfile histogram.
TProfile * BookTProfile (const std::string &name, const std::string &labels, int nBinsX, float *binsX, bool prependDir=true)
 Book a TProfile histogram with variable binning in x-axis.
TProfile * BookTProfileRangeY (const std::string &name, const std::string &labels, int nBinsX, double *binsX, double startY, double endY, bool prependDir=true)
 Book a TProfile histogram with variable binning in x-axis and limits in y-values.
TProfile2D * BookTProfile2D (const std::string &name, const std::string &labels, const int nBinsX, const double xlo, const double xhi, const int nBinsY, const double ylo, const double yhi, bool prependDir=true, bool useRMS=false)
 Book a TProfile 2D histogram with variable binning in x-axis and limits in y-values.
TProfile2D * BookTProfile2D (const std::string &name, const std::string &labels, const int nBinsX, double *binsX, const int nBinsY, double *binsY, bool prependDir=true, bool useRMS=false)
 Book a TProfile 2D histogram with variable binning in x-axis and limits in y-values.
TEfficiency * BookTEfficiency (const std::string &name, const std::string &labels, const int nBinsX, const float xlo, const float xhi, const bool prependDir=true)
 Book a (1-D) TEfficiency histogram.
TEfficiency * BookTEfficiency (const std::string &name, const std::string &labels, const int nBinsX, const float xlo, const float xhi, const int nBinsy, const float ylo, const float yhi, const bool prependDir=true)
 Book a (2-D) TEfficiency histogram.

Public Attributes

TH1 * Res_pT {nullptr}
TH1 * Res_eta {nullptr}
TH1 * Res_phi {nullptr}
TH2 * Res_pT_vs_pT {nullptr}
TH2 * Res_eta_vs_pT {nullptr}
TH2 * Res_phi_vs_pT {nullptr}
std::vector< TH2 * > Res_pT_vs_eta
std::vector< TH2 * > Res_pT_vs_phi
std::vector< TH2 * > Res_eta_vs_eta
std::vector< TH2 * > Res_phi_vs_phi

Protected Attributes

std::vector< PlotBase * > m_vSubNodes
std::vector< HistDatam_vBookedHistograms
std::vector< TreeDatam_vBookedTrees
std::vector< EfficiencyDatam_vBookedEfficiencies
std::string m_sDirectory
int m_iDetailLevel

Private Member Functions

virtual void initializePlots ()
virtual void finalizePlots ()

Static Private Member Functions

static std::string constructPrefix (std::string dir, bool prependDir)

Private Attributes

std::string m_sType
std::vector< std::string > m_pt_slices
std::vector< double > m_etaBins

Detailed Description

Definition at line 15 of file ResoTriggerMuonPlots.h.

Constructor & Destructor Documentation

◆ ResoTriggerMuonPlots()

ResoTriggerMuonPlots::ResoTriggerMuonPlots ( PlotBase * pParent,
const std::string & sDir,
std::string sType = "" )

Definition at line 13 of file ResoTriggerMuonPlots.cxx.

13 :
14 PlotBase(pParent, sDir),
15 m_sType(std::move(sType)),
16 m_pt_slices {"0", "25", "55", "100"},
17 m_etaBins {-3., -2.5, -2.4, -1.918, -1.623, -1.348, -1.2329, -1.1479, -1.05, -0.908, -0.791,
18 -0.652, -0.476, -0.324, -0.132, 0, 0.132, 0.324, 0.476, 0.652, 0.791, 0.908,
19 1.05, 1.1479, 1.2329, 1.348, 1.623, 1.918, 2.4, 2.5, 3.}
20{
21 Res_pT = Book1D("Res" + m_sType + "_pT", "Res" + m_sType + "_pT;(1/pT-1/RECOpT)/(1/RECOpT);Entries", 200, -0.25, 0.25);
22 Res_eta = Book1D("Res" + m_sType + "_eta", "Res" + m_sType + "_eta;eta-RECOeta;Entries", 200, -0.02, 0.02);
23 Res_phi = Book1D("Res" + m_sType + "_phi", "Res" + m_sType + "_phi;phi-RECOphi;Entries", 200, -0.005, 0.005);
24
25 Res_pT_vs_pT = Book2D("Res" + m_sType + "_pT_vs_pT", "Res" + m_sType + "_pT vs pT;pT [GeV];(1/pT-1/RECOpT)/(1/RECOpT)", 150, 0.,
26 150., 100, -0.25, 0.25);
28 Book2D("Res" + m_sType + "_eta_vs_pT", "Res" + m_sType + "_eta vs pT;pT [GeV];eta-RECOeta", 150, 0., 150., 200, -0.02, 0.02);
30 Book2D("Res" + m_sType + "_phi_vs_pT", "Res" + m_sType + "_phi vs pT;pT [GeV];phi-RECOphi", 150, 0., 150., 200, -0.005, 0.005);
31 Res_pT_vs_eta.clear();
32 Res_pT_vs_phi.clear();
33 Res_phi_vs_phi.clear();
34 Res_eta_vs_eta.clear();
35
36 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
37 Res_pT_vs_eta.push_back(Book2D("Res" + m_sType + "_pT_vs_eta_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
38 "Res" + m_sType + "_pT vs eta from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] +
39 " GeV; eta;(1/pT-1/RECOpT)/(1/RECOpT)",
40 m_etaBins.size() - 1, &m_etaBins[0], 100, -0.25, 0.25));
41 Res_pT_vs_phi.push_back(Book2D("Res" + m_sType + "_pT_vs_phi_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
42 "Res" + m_sType + "_pT vs phi from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] +
43 " GeV;phi;(1/pT-1/RECOpT)/(1/RECOpT)",
44 96, -M_PI, M_PI, 100, -0.25, 0.25));
45 Res_phi_vs_phi.push_back(Book2D("Res" + m_sType + "_phi_vs_phi_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
46 "Res" + m_sType + "_phi vs phi from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] +
47 " GeV;phi;(1/pT-1/RECOpT)/(1/RECOpT)",
48 96, -M_PI, M_PI, 100, -0.25, 0.25));
49 Res_eta_vs_eta.push_back(
50 Book2D("Res" + m_sType + "_eta_vs_eta_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
51 "Res" + m_sType + "_eta vs eta from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] + " GeV;eta;eta-RECOeta",
52 m_etaBins.size() - 1, &m_etaBins[0], 100, -0.25, 0.25));
53 }
54}
#define M_PI
TH1D * Book1D(const std::string &name, const std::string &labels, int nBins, float start, float end, bool prependDir=true)
Book a TH1D histogram.
Definition PlotBase.cxx:94
PlotBase(PlotBase *parent, const std::string &sDir)
Definition PlotBase.cxx:29
TH2F * Book2D(const std::string &name, const std::string &labels, int nBinsX, float startX, float endX, int nBinsY, float startY, float endY, bool prependDir=true)
Book a TH2F histogram.
Definition PlotBase.cxx:123
std::vector< double > m_etaBins
std::vector< TH2 * > Res_eta_vs_eta
std::vector< TH2 * > Res_pT_vs_eta
std::vector< TH2 * > Res_pT_vs_phi
std::vector< TH2 * > Res_phi_vs_phi
std::vector< std::string > m_pt_slices

Member Function Documentation

◆ Book1D() [1/2]

TH1D * PlotBase::Book1D ( const std::string & name,
const std::string & labels,
int nBins,
float start,
float end,
bool prependDir = true )
inherited

Book a TH1D histogram.

Definition at line 94 of file PlotBase.cxx.

95 {
96 std::string prefix = constructPrefix(m_sDirectory, prependDir);
97 Bool_t oldstat = TH1::AddDirectoryStatus();
98 TH1::AddDirectory(false);
99 TH1D *hist = new TH1D((prefix + name).c_str(), labels.c_str(), nBins, start, end);
100 TH1::AddDirectory(oldstat);
101
102 hist->Sumw2();
103 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
104 return hist;
105}
static std::string constructPrefix(std::string dir, bool prependDir)
Definition PlotBase.cxx:293
std::vector< HistData > m_vBookedHistograms
Definition PlotBase.h:97
std::string m_sDirectory
Definition PlotBase.h:100

◆ Book1D() [2/2]

TH1D * PlotBase::Book1D ( const std::string & name,
TH1 * refHist,
const std::string & labels,
bool prependDir = true )
inherited

Book a TH1D histogram using refHist as reference for number of bins and axis range.

Definition at line 108 of file PlotBase.cxx.

108 {
109 std::string prefix = constructPrefix(m_sDirectory, prependDir);
110 Bool_t oldstat = TH1::AddDirectoryStatus();
111 TH1::AddDirectory(false);
112 TH1D *hist = new TH1D((prefix + name).c_str(), labels.c_str(), refHist->GetNbinsX(),
113 refHist->GetXaxis()->GetXbins()->GetArray());
114 hist->Sumw2();
115 TH1::AddDirectory(oldstat);
116
117
118 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
119 return hist;
120}

◆ Book2D() [1/3]

TH2F * PlotBase::Book2D ( const std::string & name,
const std::string & labels,
int nBinsX,
Double_t * binsX,
int nBinsY,
Double_t startY,
Double_t endY,
bool prependDir = true )
inherited

Book a TH2F histogram with variable x axis binning.

Definition at line 144 of file PlotBase.cxx.

145 {
146 std::string prefix = constructPrefix(m_sDirectory, prependDir);
147 Bool_t oldstat = TH2::AddDirectoryStatus();
148 TH2::AddDirectory(false);
149 TH2F *hist = new TH2F((prefix + name).c_str(), labels.c_str(), nBinsX, binsX, nBinsY, startY, endY);
150 hist->Sumw2();
151 TH2::AddDirectory(oldstat);
152 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
153 return hist;
154}
TH2F(name, title, nxbins, bins_par2, bins_par3, bins_par4, bins_par5=None, bins_par6=None, path='', **kwargs)

◆ Book2D() [2/3]

TH2F * PlotBase::Book2D ( const std::string & name,
const std::string & labels,
int nBinsX,
float startX,
float endX,
int nBinsY,
float startY,
float endY,
bool prependDir = true )
inherited

Book a TH2F histogram.

Definition at line 123 of file PlotBase.cxx.

124 {
125 std::string prefix = constructPrefix(m_sDirectory, prependDir);
126 Bool_t oldstat = TH2::AddDirectoryStatus();
127 TH2::AddDirectory(false);
128 TH2F *hist = new TH2F((prefix + name).c_str(), labels.c_str(), nBinsX, startX, endX, nBinsY, startY, endY);
129 hist->Sumw2();
130 TH2::AddDirectory(oldstat);
131
132
133 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
134 return hist;
135}

◆ Book2D() [3/3]

TH2F * PlotBase::Book2D ( const std::string & name,
TH2 * refHist,
const std::string & labels,
bool prependDir = true )
inherited

Book a TH2D histogram using refHist as reference for number of bins and axis range.

Definition at line 138 of file PlotBase.cxx.

138 {
139 return Book2D(name, labels, refHist->GetNbinsX(), refHist->GetXaxis()->GetXmin(), refHist->GetXaxis()->GetXmax(),
140 refHist->GetNbinsY(), refHist->GetYaxis()->GetXmin(), refHist->GetYaxis()->GetXmax(), prependDir);
141}

◆ Book3D() [1/2]

TH3F * PlotBase::Book3D ( const std::string & name,
const std::string & labels,
int nBinsX,
float startX,
float endX,
int nBinsY,
float startY,
float endY,
int nBinsZ,
float startZ,
float endZ,
bool prependDir = true )
inherited

Book a TH3F histogram.

Definition at line 157 of file PlotBase.cxx.

158 {
159 std::string prefix = constructPrefix(m_sDirectory, prependDir);
160 Bool_t oldstat = TH3::AddDirectoryStatus();
161 TH3::AddDirectory(false);
162 TH3F *hist = new TH3F((prefix + name).c_str(),
163 labels.c_str(), nBinsX, startX, endX, nBinsY, startY, endY, nBinsZ, startZ, endZ);
164 hist->Sumw2();
165 TH3::AddDirectory(oldstat);
166 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
167 return hist;
168}

◆ Book3D() [2/2]

TH3F * PlotBase::Book3D ( const std::string & name,
TH3 * refHist,
const std::string & labels,
bool prependDir = true )
inherited

Book a TH3F histogram using refHist as reference for number of bins and axis range.

Definition at line 171 of file PlotBase.cxx.

171 {
172 std::string prefix = constructPrefix(m_sDirectory, prependDir);
173 Bool_t oldstat = TH3::AddDirectoryStatus();
174 TH3::AddDirectory(false);
175 TH3F *hist = new TH3F((prefix + name).c_str(), labels.c_str(), refHist->GetNbinsX(),
176 refHist->GetXaxis()->GetXbins()->GetArray(), refHist->GetNbinsY(),
177 refHist->GetYaxis()->GetXbins()->GetArray(), refHist->GetNbinsZ(),
178 refHist->GetZaxis()->GetXbins()->GetArray());
179 TH3::AddDirectory(oldstat);
180
181 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
182 return hist;
183}

◆ BookTEfficiency() [1/2]

TEfficiency * PlotBase::BookTEfficiency ( const std::string & name,
const std::string & labels,
const int nBinsX,
const float xlo,
const float xhi,
const bool prependDir = true )
inherited

Book a (1-D) TEfficiency histogram.

Definition at line 257 of file PlotBase.cxx.

257 {
258 std::string prefix = constructPrefix(m_sDirectory, prependDir);
259 //Bool_t oldstat = TEfficiency::AddDirectoryStatus();
260 TEfficiency *hist = new TEfficiency((prefix + name).c_str(), labels.c_str(), nBinsX, xlo, xhi);
261 //hist->SetAutoSave(0);
262 //hist->SetAtoFlush(0);
263 hist->SetDirectory(nullptr);
264 m_vBookedEfficiencies.emplace_back(hist, m_sDirectory);
265 //TEfficiency::AddDirectory(oldstat);
266 return hist;
267}
std::vector< EfficiencyData > m_vBookedEfficiencies
Definition PlotBase.h:99

◆ BookTEfficiency() [2/2]

TEfficiency * PlotBase::BookTEfficiency ( const std::string & name,
const std::string & labels,
const int nBinsX,
const float xlo,
const float xhi,
const int nBinsy,
const float ylo,
const float yhi,
const bool prependDir = true )
inherited

Book a (2-D) TEfficiency histogram.

Definition at line 270 of file PlotBase.cxx.

270 {
271 std::string prefix = constructPrefix(m_sDirectory, prependDir);
272
273 TEfficiency *hist = new TEfficiency((prefix + name).c_str(), labels.c_str(), nBinsX, xlo, xhi, nBinsY, ylo, yhi);
274 hist->SetDirectory(nullptr);
275 m_vBookedEfficiencies.emplace_back(hist, m_sDirectory);
276
277 return hist;
278}

◆ BookTProfile() [1/2]

TProfile * PlotBase::BookTProfile ( const std::string & name,
const std::string & labels,
int nBinsX,
float * binsX,
bool prependDir = true )
inherited

Book a TProfile histogram with variable binning in x-axis.

Definition at line 204 of file PlotBase.cxx.

204 {
205 std::string prefix = constructPrefix(m_sDirectory, prependDir);
206 TProfile *hist(nullptr);
207 Bool_t oldstat = TProfile::AddDirectoryStatus();
208 TProfile::AddDirectory(false);
209
210 hist = new TProfile((prefix + name).c_str(), labels.c_str(), nBinsX, binsX);
211 TProfile::AddDirectory(oldstat);
212 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
213 return hist;
214}

◆ BookTProfile() [2/2]

TProfile * PlotBase::BookTProfile ( const std::string & name,
const std::string & labels,
int nBinsX,
float startX,
float endX,
float startY = -1,
float endY = -1,
bool prependDir = true,
bool useRMS = false )
inherited

Book a TProfile histogram.

Definition at line 186 of file PlotBase.cxx.

187 {
188 std::string prefix = constructPrefix(m_sDirectory, prependDir);
189 TProfile *hist(nullptr);
190 Bool_t oldstat = TProfile::AddDirectoryStatus();
191 TProfile::AddDirectory(false);
192 std::string opt = useRMS ? "S" : "";
193 if ((startY == -1) and (endY == -1)) {
194 hist = new TProfile((prefix + name).c_str(), labels.c_str(), nBinsX, startX, endX, opt.c_str());
195 } else {
196 hist = new TProfile((prefix + name).c_str(), labels.c_str(), nBinsX, startX, endX, startY, endY, opt.c_str());
197 }
198 TProfile::AddDirectory(oldstat);
199 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
200 return hist;
201}

◆ BookTProfile2D() [1/2]

TProfile2D * PlotBase::BookTProfile2D ( const std::string & name,
const std::string & labels,
const int nBinsX,
const double xlo,
const double xhi,
const int nBinsY,
const double ylo,
const double yhi,
bool prependDir = true,
bool useRMS = false )
inherited

Book a TProfile 2D histogram with variable binning in x-axis and limits in y-values.

Definition at line 231 of file PlotBase.cxx.

233 {
234 std::string prefix = constructPrefix(m_sDirectory, prependDir);
235 Bool_t oldstat = TProfile2D::AddDirectoryStatus();
236 TProfile2D::AddDirectory(false);
237 std::string opt = useRMS ? "S" : "";
238 TProfile2D *hist = new TProfile2D((prefix + name).c_str(), labels.c_str(), nBinsX, xlo, xhi, nBinsY, ylo, yhi, opt.c_str());
239 TProfile2D::AddDirectory(oldstat);
240 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
241 return hist;
242}

◆ BookTProfile2D() [2/2]

TProfile2D * PlotBase::BookTProfile2D ( const std::string & name,
const std::string & labels,
const int nBinsX,
double * binsX,
const int nBinsY,
double * binsY,
bool prependDir = true,
bool useRMS = false )
inherited

Book a TProfile 2D histogram with variable binning in x-axis and limits in y-values.

Definition at line 245 of file PlotBase.cxx.

245 {
246 std::string prefix = constructPrefix(m_sDirectory, prependDir);
247 Bool_t oldstat = TProfile2D::AddDirectoryStatus();
248 TProfile2D::AddDirectory(false);
249 std::string opt = useRMS ? "S" : "";
250 TProfile2D *hist = new TProfile2D((prefix + name).c_str(), labels.c_str(), nBinsX, binsX, nBinsY, binsY, opt.c_str());
251 TProfile2D::AddDirectory(oldstat);
252 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
253 return hist;
254}

◆ BookTProfileRangeY()

TProfile * PlotBase::BookTProfileRangeY ( const std::string & name,
const std::string & labels,
int nBinsX,
double * binsX,
double startY,
double endY,
bool prependDir = true )
inherited

Book a TProfile histogram with variable binning in x-axis and limits in y-values.

Definition at line 217 of file PlotBase.cxx.

218 {
219 std::string prefix = constructPrefix(m_sDirectory, prependDir);
220 TProfile *hist(nullptr);
221 Bool_t oldstat = TProfile::AddDirectoryStatus();
222 TProfile::AddDirectory(false);
223
224 hist = new TProfile((prefix + name).c_str(), labels.c_str(), (Int_t) nBinsX, binsX, startY, endY);
225 TProfile::AddDirectory(oldstat);
226 m_vBookedHistograms.emplace_back(hist, m_sDirectory);
227 return hist;
228}

◆ BookTree()

TTree * PlotBase::BookTree ( const std::string & name,
bool prependDir = true )
inherited

Book a TTree.

Definition at line 281 of file PlotBase.cxx.

281 {
282 std::string prefix = constructPrefix(m_sDirectory, prependDir);
283 TTree *tree = new TTree((prefix + name).c_str(), "");
284
285 tree->SetAutoSave(0);
286 tree->SetAutoFlush(0);
287 tree->SetDirectory(nullptr);
288 m_vBookedTrees.emplace_back(tree, m_sDirectory);
289 return tree;
290}
std::vector< TreeData > m_vBookedTrees
Definition PlotBase.h:98
TChain * tree

◆ constructPrefix()

std::string PlotBase::constructPrefix ( std::string dir,
bool prependDir )
staticprivateinherited

Definition at line 293 of file PlotBase.cxx.

293 {
294 if (!prependDir) {
295 return "";
296 }
297 std::replace(dir.begin(), dir.end(), '/', '_');
298 return dir;
299}

◆ fill() [1/3]

void ResoTriggerMuonPlots::fill ( const xAOD::L2CombinedMuon & L2CBmu,
const xAOD::Muon & Recomu )

Definition at line 116 of file ResoTriggerMuonPlots.cxx.

116 {
117 float pt = 0.001 * Recomu.pt();
118 float eta = Recomu.eta();
119 float phi = Recomu.phi();
120 float respt = (1. / (L2CBmu.pt()) - 1. / Recomu.pt()) / (1. / Recomu.pt());
121 ;
122
123 const float d_phi = deltaPhi(L2CBmu.phi() , phi);
124 const float d_eta = L2CBmu.eta() - eta;
125 Res_pT->Fill(respt);
126 Res_eta->Fill(d_eta);
127 Res_phi->Fill(d_phi);
128
129 Res_pT_vs_pT->Fill(pt, respt);
130 Res_eta_vs_pT->Fill(pt, d_eta);
131 Res_phi_vs_pT->Fill(pt, d_phi);
132
133 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
134 float pt_min = atof(m_pt_slices[i].c_str());
135 float pt_max = atof(m_pt_slices[i + 1].c_str());
136 if (pt >= pt_min && pt < pt_max) {
137 Res_phi_vs_phi[i]->Fill(phi, d_phi);
138 Res_eta_vs_eta[i]->Fill(eta, d_eta);
139 Res_pT_vs_eta[i]->Fill(eta, respt);
140 Res_pT_vs_phi[i]->Fill(phi, respt);
141 }
142 }
143}
Scalar eta() const
pseudorapidity method
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
Scalar phi() const
phi method
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
double atof(std::string_view str)
Converts a string into a double / float.

◆ fill() [2/3]

void ResoTriggerMuonPlots::fill ( const xAOD::L2StandAloneMuon & L2SAmu,
const xAOD::Muon & Recomu )

Definition at line 86 of file ResoTriggerMuonPlots.cxx.

86 {
87 float pt = 0.001 * Recomu.pt();
88 float eta = Recomu.eta();
89 float phi = Recomu.phi();
90 float respt = 0;
91 if (L2SAmu.pt() > 0.) respt = (1. / (1000. * L2SAmu.pt()) - 1. / Recomu.pt()) / (1. / Recomu.pt());
92 if (L2SAmu.pt() < 0.) respt = (1. / (-1000. * L2SAmu.pt()) - 1. / Recomu.pt()) / (1. / Recomu.pt());
93
94 const float d_eta = L2SAmu.eta() - eta;
95 const float d_phi = deltaPhi(L2SAmu.phi(), phi);
96 Res_pT->Fill(respt);
97 Res_eta->Fill(d_eta);
98 Res_phi->Fill(d_phi);
99
100 Res_pT_vs_pT->Fill(pt, respt);
101 Res_eta_vs_pT->Fill(pt, d_eta);
102 Res_phi_vs_pT->Fill(pt, d_phi);
103
104 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
105 float pt_min = atof(m_pt_slices[i].c_str());
106 float pt_max = atof(m_pt_slices[i + 1].c_str());
107 if (pt >= pt_min && pt < pt_max) {
108 Res_phi_vs_phi[i]->Fill(phi, d_phi);
109 Res_eta_vs_eta[i]->Fill(eta, d_eta);
110 Res_pT_vs_eta[i]->Fill(eta, respt);
111 Res_pT_vs_phi[i]->Fill(phi, respt);
112 }
113 }
114}
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double eta() const
The pseudorapidity ( ) of the particle.

◆ fill() [3/3]

void ResoTriggerMuonPlots::fill ( const xAOD::Muon & Trigmu,
const xAOD::Muon & Recomu )

Definition at line 56 of file ResoTriggerMuonPlots.cxx.

56 {
57 // fill plots Res_pt,eta,phi
58 float pt = 0.001 * Recomu.pt();
59 float eta = Recomu.eta();
60 float phi = Recomu.phi();
61 // float respt = (Trigmu.pt() - Recomu.pt())/Recomu.pt();
62 float respt = (1. / Trigmu.pt() - 1. / Recomu.pt()) / (1. / Recomu.pt());
63
64 const float d_phi = deltaPhi(Trigmu.phi(), phi);
65 const float d_eta = Trigmu.eta() - eta;
66 Res_pT->Fill(respt);
67 Res_eta->Fill(d_eta);
68 Res_phi->Fill(d_phi);
69
70 Res_pT_vs_pT->Fill(pt, respt);
71 Res_eta_vs_pT->Fill(pt, d_eta);
72 Res_phi_vs_pT->Fill(pt, d_phi);
73
74 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
75 float pt_min = atof(m_pt_slices[i].c_str());
76 float pt_max = atof(m_pt_slices[i + 1].c_str());
77 if (pt >= pt_min && pt < pt_max) {
78 Res_pT_vs_eta[i]->Fill(eta, respt);
79 Res_pT_vs_phi[i]->Fill(phi, respt);
80 Res_phi_vs_phi[i]->Fill(phi, d_phi);
81 Res_eta_vs_eta[i]->Fill(eta, d_eta);
82 }
83 }
84}

◆ finalize()

void PlotBase::finalize ( )
inherited

Definition at line 47 of file PlotBase.cxx.

47 {
48 for (auto *subNode: m_vSubNodes) {
49 subNode->finalize();
50 }
52}
std::vector< PlotBase * > m_vSubNodes
Definition PlotBase.h:96
virtual void finalizePlots()
Definition PlotBase.h:92

◆ finalizePlots()

◆ getDirectory()

const std::string & PlotBase::getDirectory ( )
inlineinherited

Definition at line 88 of file PlotBase.h.

88{return m_sDirectory;}

◆ initialize()

void PlotBase::initialize ( )
inherited

Definition at line 39 of file PlotBase.cxx.

39 {
40 for (auto *subNode: m_vSubNodes) {
41 subNode->initialize();
42 }
44}
virtual void initializePlots()
Definition PlotBase.h:91

◆ initializePlots()

virtual void PlotBase::initializePlots ( )
inlineprivatevirtualinherited

Reimplemented in DiTau::CorePlots, DiTau::ResolutionPlots, Egamma::ClusMomentumPlots, Egamma::ElectronFrwdPlots, Egamma::ElectronPlots, Egamma::IsolationPlots, Egamma::KinematicsPlots, Egamma::LRTElectronPlots, Egamma::PhotonAmbPlots, Egamma::PhotonCnvPlots, Egamma::PhotonConversionPlots, Egamma::PhotonPlots, Egamma::ShowerShapesPlots, Egamma::TrackPlots, ElectronValidationPlots, IDTPM::DuplicateRatePlots, IDTPM::EfficiencyPlots, IDTPM::FakeRatePlots, IDTPM::HitsOnTracksPlots, IDTPM::NtracksPlots, IDTPM::OfflineElectronPlots, IDTPM::ResolutionPlots, IDTPM::SummaryPlots, IDTPM::TrackParametersPlots, IDTPM::VertexParametersPlots, InDetBasicPlot, InDetPerfNtuple, InDetPerfPlot_Duplicate, InDetPerfPlot_Efficiency, InDetPerfPlot_FakeRate, InDetPerfPlot_HitEfficiency, InDetPerfPlot_HitResidual, InDetPerfPlot_Hits, InDetPerfPlot_nTracks, InDetPerfPlot_Resolution, InDetPerfPlot_TrackParameters, InDetPerfPlot_TrkInJet, InDetPerfPlot_TRTExtension, InDetPerfPlot_Vertex, InDetPerfPlot_VertexTruthMatching, InDetPerfPlot_VerticesVsMu, JetTagDQA::BTaggingValidationPlots, LRTElectronValidationPlots, Muon::BetaPlots, Muon::ChargeDepParamPlots, Muon::HitFracTypePlots, Muon::IsoCorrPlots, Muon::IsoPlots, Muon::MomentumPullPlots, Muon::MomentumTruthPullPlots, Muon::MuonParamElossPlots, Muon::MuonParamPlots, Muon::MuonTree, Muon::RecoInfoPlots, Muon::SlowMuonParamPlots, PFO::ClusterMomentPlots, PFO::ClusterPlots, PFO::FlowElement_LinkerPlots, PFO::LeptonFELinkerPlots, PFO::PFOAlgPropertyPlots, PFO::PFOAttributePlots, PFO::PFOCalibHitClusterMomentPlots, PFO::PFOClusterMomentPlots, PFO::PFOPlots, PFO::PFOPVMatchedPlots, PhotonValidationPlots, PhysVal::BTagPlots, PhysVal::EventInfoPlots, PhysVal::KinematicsPlots, PhysVal::METPlots, PhysVal::TrkAndVtxPlots, RecoLumiPlots, RecoMuonIDTrackPlots, RecoMuonPlots, RecoMuonSegmentPlots, RecoMuonTrackPlots, RecoPhysPlots, RecoVertexPlots, Tau::CorePlots, Tau::DecayModeMigration, Tau::EfficiencyPlots, Tau::EVetoPlots, Tau::GeneralTauPlots, Tau::ResolutionPlots, Tau::TauIDVariablesPlots, Tau::TauKinematicPlots, Tau::TauParticleFlowPlots, TCCPlots, Trk::DefParamPullPlots, Trk::EfficiencyPlots, Trk::ExtrLayerPlots, Trk::ExtrRegionPlots, Trk::HitResidualPlots, Trk::HitTypePlots, Trk::IDHitPlots, Trk::ImpactPlots, Trk::ParamPlots, Trk::RecoInfoPlots, Trk::ResolutionPlots, Trk::TruthInfoPlots, ZeeValidation::FWDZeePlots, ZeeValidation::ReconElectronsPlots, ZeeValidation::TrueElectronsPlots, ZeeValidation::TrueFwdElectronsPlots, and ZeeValidation::ZeePlots.

Definition at line 91 of file PlotBase.h.

91{;}

◆ RegisterSubPlot()

void PlotBase::RegisterSubPlot ( PlotBase * pPlotBase)
inlineinherited

Definition at line 41 of file PlotBase.h.

41{m_vSubNodes.push_back(pPlotBase);}

◆ retrieveBookedEfficiencies()

std::vector< EfficiencyData > PlotBase::retrieveBookedEfficiencies ( )
inherited

Retrieve all booked efficiency objects.

Definition at line 83 of file PlotBase.cxx.

83 {
84 std::vector<EfficiencyData> vBookedEfficiencies = m_vBookedEfficiencies;
85 for (const auto &subNode: m_vSubNodes) {
86 std::vector<EfficiencyData> subNodeHists = subNode->retrieveBookedEfficiencies();
87 vBookedEfficiencies.insert(vBookedEfficiencies.end(), subNodeHists.begin(), subNodeHists.end());
88 }
89 return vBookedEfficiencies;
90}

◆ retrieveBookedHistograms()

std::vector< HistData > PlotBase::retrieveBookedHistograms ( )
inherited

Retrieve all booked histograms.

Definition at line 63 of file PlotBase.cxx.

63 {
64 std::vector<HistData> vBookedHistograms = m_vBookedHistograms;
65 for (const auto &subNode: m_vSubNodes) {
66 std::vector<HistData> subNodeHists = subNode->retrieveBookedHistograms();
67 vBookedHistograms.insert(vBookedHistograms.end(), subNodeHists.begin(), subNodeHists.end());
68 }
69 return vBookedHistograms;
70}

◆ retrieveBookedTrees()

std::vector< TreeData > PlotBase::retrieveBookedTrees ( )
inherited

Retrieve all booked trees.

Definition at line 73 of file PlotBase.cxx.

73 {
74 std::vector<TreeData> vBookedTrees = m_vBookedTrees;
75 for (auto *subNode: m_vSubNodes) {
76 std::vector<TreeData> subNodeTrees = subNode->retrieveBookedTrees();
77 vBookedTrees.insert(vBookedTrees.end(), subNodeTrees.begin(), subNodeTrees.end());
78 }
79 return vBookedTrees;
80}

◆ setDetailLevel()

void PlotBase::setDetailLevel ( int iDetailLevel)
inherited

Definition at line 55 of file PlotBase.cxx.

55 {
56 for (auto *subNode: m_vSubNodes) {
57 subNode->setDetailLevel(iDetailLevel);
58 }
59 m_iDetailLevel = iDetailLevel;
60}
int m_iDetailLevel
Definition PlotBase.h:101

Member Data Documentation

◆ m_etaBins

std::vector<double> ResoTriggerMuonPlots::m_etaBins
private

Definition at line 36 of file ResoTriggerMuonPlots.h.

◆ m_iDetailLevel

int PlotBase::m_iDetailLevel
protectedinherited

Definition at line 101 of file PlotBase.h.

◆ m_pt_slices

std::vector<std::string> ResoTriggerMuonPlots::m_pt_slices
private

Definition at line 35 of file ResoTriggerMuonPlots.h.

◆ m_sDirectory

std::string PlotBase::m_sDirectory
protectedinherited

Definition at line 100 of file PlotBase.h.

◆ m_sType

std::string ResoTriggerMuonPlots::m_sType
private

Definition at line 34 of file ResoTriggerMuonPlots.h.

◆ m_vBookedEfficiencies

std::vector<EfficiencyData> PlotBase::m_vBookedEfficiencies
protectedinherited

Definition at line 99 of file PlotBase.h.

◆ m_vBookedHistograms

std::vector<HistData> PlotBase::m_vBookedHistograms
protectedinherited

Definition at line 97 of file PlotBase.h.

◆ m_vBookedTrees

std::vector<TreeData> PlotBase::m_vBookedTrees
protectedinherited

Definition at line 98 of file PlotBase.h.

◆ m_vSubNodes

std::vector<PlotBase*> PlotBase::m_vSubNodes
protectedinherited

Definition at line 96 of file PlotBase.h.

◆ Res_eta

TH1* ResoTriggerMuonPlots::Res_eta {nullptr}

Definition at line 23 of file ResoTriggerMuonPlots.h.

23{nullptr};

◆ Res_eta_vs_eta

std::vector<TH2*> ResoTriggerMuonPlots::Res_eta_vs_eta

Definition at line 30 of file ResoTriggerMuonPlots.h.

◆ Res_eta_vs_pT

TH2* ResoTriggerMuonPlots::Res_eta_vs_pT {nullptr}

Definition at line 26 of file ResoTriggerMuonPlots.h.

26{nullptr};

◆ Res_phi

TH1* ResoTriggerMuonPlots::Res_phi {nullptr}

Definition at line 24 of file ResoTriggerMuonPlots.h.

24{nullptr};

◆ Res_phi_vs_phi

std::vector<TH2*> ResoTriggerMuonPlots::Res_phi_vs_phi

Definition at line 31 of file ResoTriggerMuonPlots.h.

◆ Res_phi_vs_pT

TH2* ResoTriggerMuonPlots::Res_phi_vs_pT {nullptr}

Definition at line 27 of file ResoTriggerMuonPlots.h.

27{nullptr};

◆ Res_pT

TH1* ResoTriggerMuonPlots::Res_pT {nullptr}

Definition at line 22 of file ResoTriggerMuonPlots.h.

22{nullptr};

◆ Res_pT_vs_eta

std::vector<TH2*> ResoTriggerMuonPlots::Res_pT_vs_eta

Definition at line 28 of file ResoTriggerMuonPlots.h.

◆ Res_pT_vs_phi

std::vector<TH2*> ResoTriggerMuonPlots::Res_pT_vs_phi

Definition at line 29 of file ResoTriggerMuonPlots.h.

◆ Res_pT_vs_pT

TH2* ResoTriggerMuonPlots::Res_pT_vs_pT {nullptr}

Definition at line 25 of file ResoTriggerMuonPlots.h.

25{nullptr};

The documentation for this class was generated from the following files: