Analysis Software
Documentation for
sPHENIX
simulation software
Home page
Related Pages
Modules
Namespaces
Classes
Files
Examples
External Links
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
FakeRatePlotTool.cpp
Go to the documentation of this file.
Or view
the newest version in sPHENIX GitHub for file FakeRatePlotTool.cpp
1
// This file is part of the Acts project.
2
//
3
// Copyright (C) 2019 CERN for the benefit of the Acts project
4
//
5
// This Source Code Form is subject to the terms of the Mozilla Public
6
// License, v. 2.0. If a copy of the MPL was not distributed with this
7
// file, You can obtain one at http://mozilla.org/MPL/2.0/.
8
9
#include "
ActsExamples/Validation/FakeRatePlotTool.hpp
"
10
11
#include "
Acts/Definitions/Algebra.hpp
"
12
#include "
Acts/Utilities/VectorHelpers.hpp
"
13
#include "
ActsFatras/EventData/Particle.hpp
"
14
15
#include <TEfficiency.h>
16
#include <TH2.h>
17
18
using
Acts::VectorHelpers::eta
;
19
using
Acts::VectorHelpers::perp
;
20
using
Acts::VectorHelpers::phi
;
21
using
Acts::VectorHelpers::theta
;
22
23
ActsExamples::FakeRatePlotTool::FakeRatePlotTool
(
24
const
ActsExamples::FakeRatePlotTool::Config
&
cfg
,
Acts::Logging::Level
lvl)
25
:
m_cfg
(cfg),
m_logger
(Acts::
getDefaultLogger
(
"FakeRatePlotTool"
, lvl)) {}
26
27
void
ActsExamples::FakeRatePlotTool::book
(
28
FakeRatePlotTool::FakeRatePlotCache
& fakeRatePlotCache)
const
{
29
PlotHelpers::Binning
bPt =
m_cfg
.varBinning.at(
"Pt"
);
30
PlotHelpers::Binning
bEta =
m_cfg
.varBinning.at(
"Eta"
);
31
PlotHelpers::Binning
bPhi =
m_cfg
.varBinning.at(
"Phi"
);
32
PlotHelpers::Binning
bNum =
m_cfg
.varBinning.at(
"Num"
);
33
ACTS_DEBUG
(
"Initialize the histograms for fake rate plots"
);
34
35
// number of reco tracks vs pT scatter plots
36
fakeRatePlotCache.
nReco_vs_pT
=
PlotHelpers::bookHisto
(
37
"nRecoTracks_vs_pT"
,
"Number of reconstructed track candidates"
, bPt,
38
bNum);
39
// number of truth-matched tracks vs pT scatter plots
40
fakeRatePlotCache.
nTruthMatched_vs_pT
=
PlotHelpers::bookHisto
(
41
"nTruthMatchedTracks_vs_pT"
,
"Number of truth-matched track candidates"
,
42
bPt, bNum);
43
// number of fake tracks vs pT scatter plots
44
fakeRatePlotCache.
nFake_vs_pT
=
PlotHelpers::bookHisto
(
45
"nFakeTracks_vs_pT"
,
"Number of fake track candidates"
, bPt, bNum);
46
47
// number of reco tracks vs eta scatter plots
48
fakeRatePlotCache.
nReco_vs_eta
=
PlotHelpers::bookHisto
(
49
"nRecoTracks_vs_eta"
,
"Number of reconstructed track candidates"
, bEta,
50
bNum);
51
// number of truth-matched tracks vs eta scatter plots
52
fakeRatePlotCache.
nTruthMatched_vs_eta
=
PlotHelpers::bookHisto
(
53
"nTruthMatchedTracks_vs_eta"
,
"Number of truth-matched track candidates"
,
54
bEta, bNum);
55
// number of fake tracks vs eta scatter plots
56
fakeRatePlotCache.
nFake_vs_eta
=
PlotHelpers::bookHisto
(
57
"nFakeTracks_vs_eta"
,
"Number of fake track candidates"
, bEta, bNum);
58
59
// fake rate vs pT
60
fakeRatePlotCache.
fakeRate_vs_pT
=
PlotHelpers::bookEff
(
61
"fakerate_vs_pT"
,
"Tracking fake rate;pT [GeV/c];Fake rate"
, bPt);
62
// fake rate vs eta
63
fakeRatePlotCache.
fakeRate_vs_eta
=
PlotHelpers::bookEff
(
64
"fakerate_vs_eta"
,
"Tracking fake rate;#eta;Fake rate"
, bEta);
65
// fake rate vs phi
66
fakeRatePlotCache.
fakeRate_vs_phi
=
PlotHelpers::bookEff
(
67
"fakerate_vs_phi"
,
"Tracking fake rate;#phi;Fake rate"
, bPhi);
68
}
69
70
void
ActsExamples::FakeRatePlotTool::clear
(
71
FakeRatePlotCache
& fakeRatePlotCache)
const
{
72
delete
fakeRatePlotCache.
nReco_vs_pT
;
73
delete
fakeRatePlotCache.
nTruthMatched_vs_pT
;
74
delete
fakeRatePlotCache.
nFake_vs_pT
;
75
delete
fakeRatePlotCache.
nReco_vs_eta
;
76
delete
fakeRatePlotCache.
nTruthMatched_vs_eta
;
77
delete
fakeRatePlotCache.
nFake_vs_eta
;
78
delete
fakeRatePlotCache.
fakeRate_vs_pT
;
79
delete
fakeRatePlotCache.
fakeRate_vs_eta
;
80
delete
fakeRatePlotCache.
fakeRate_vs_phi
;
81
}
82
83
void
ActsExamples::FakeRatePlotTool::write
(
84
const
FakeRatePlotTool::FakeRatePlotCache
& fakeRatePlotCache)
const
{
85
ACTS_DEBUG
(
"Write the plots to output file."
);
86
fakeRatePlotCache.
nReco_vs_pT
->Write();
87
fakeRatePlotCache.
nTruthMatched_vs_pT
->Write();
88
fakeRatePlotCache.
nFake_vs_pT
->Write();
89
fakeRatePlotCache.
nReco_vs_eta
->Write();
90
fakeRatePlotCache.
nTruthMatched_vs_eta
->Write();
91
fakeRatePlotCache.
nFake_vs_eta
->Write();
92
fakeRatePlotCache.
fakeRate_vs_pT
->Write();
93
fakeRatePlotCache.
fakeRate_vs_eta
->Write();
94
fakeRatePlotCache.
fakeRate_vs_phi
->Write();
95
}
96
97
void
ActsExamples::FakeRatePlotTool::fill
(
98
FakeRatePlotTool::FakeRatePlotCache
& fakeRatePlotCache,
99
const
Acts::BoundTrackParameters
& fittedParameters,
bool
status
)
const
{
100
const
auto
momentum
= fittedParameters.
momentum
();
101
const
double
fit_phi =
phi
(
momentum
);
102
const
double
fit_eta =
eta
(
momentum
);
103
const
double
fit_pT =
perp
(
momentum
);
104
105
PlotHelpers::fillEff
(fakeRatePlotCache.
fakeRate_vs_pT
, fit_pT, status);
106
PlotHelpers::fillEff
(fakeRatePlotCache.
fakeRate_vs_eta
, fit_eta, status);
107
PlotHelpers::fillEff
(fakeRatePlotCache.
fakeRate_vs_phi
, fit_phi, status);
108
}
109
110
void
ActsExamples::FakeRatePlotTool::fill
(
111
FakeRatePlotTool::FakeRatePlotCache
& fakeRatePlotCache,
112
const
ActsFatras::Particle
& truthParticle,
size_t
nTruthMatchedTracks,
113
size_t
nFakeTracks)
const
{
114
const
auto
t_eta =
eta
(truthParticle.
direction
());
115
const
auto
t_pT = truthParticle.
transverseMomentum
();
116
117
PlotHelpers::fillHisto
(fakeRatePlotCache.
nReco_vs_pT
, t_pT,
118
nTruthMatchedTracks + nFakeTracks);
119
PlotHelpers::fillHisto
(fakeRatePlotCache.
nTruthMatched_vs_pT
, t_pT,
120
nTruthMatchedTracks);
121
PlotHelpers::fillHisto
(fakeRatePlotCache.
nFake_vs_pT
, t_pT, nFakeTracks);
122
123
PlotHelpers::fillHisto
(fakeRatePlotCache.
nReco_vs_eta
, t_eta,
124
nTruthMatchedTracks + nFakeTracks);
125
PlotHelpers::fillHisto
(fakeRatePlotCache.
nTruthMatched_vs_eta
, t_eta,
126
nTruthMatchedTracks);
127
PlotHelpers::fillHisto
(fakeRatePlotCache.
nFake_vs_eta
, t_eta, nFakeTracks);
128
}
acts
blob
sPHENIX
Examples
Framework
src
Validation
FakeRatePlotTool.cpp
Built by
Jin Huang
. updated:
Sat Feb 17 2024 22:17:38
using
1.8.2 with
sPHENIX GitHub integration