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

get resolution for electron and photons (converted / unconverted) vs E,eta More...

#include <eg_resolution.h>

Collaboration diagram for eg_resolution:

Public Member Functions

 eg_resolution (const std::string &configuration)
 constructor (initialization done there reading root files with resolution fit parameters
 ~eg_resolution ()=default
double getResolution (int particle_type, double energy, double eta, int resol_type=2) const
 get relative resolution (sigmaE/E) as a function of E (in Mev) and eta
double getResolution (const xAOD::Egamma &particle, int resol_type=2) const
 get relative resolution (sigmaE/E) for egamma particles

Static Public Attributes

static constexpr size_t samplings = 4
static constexpr size_t resolTypes = 3

Private Attributes

std::array< std::array< std::unique_ptr< TH1 >, resolTypes >, samplingsm_hSampling
std::array< std::array< std::unique_ptr< TH1 >, resolTypes >, samplingsm_hNoise
std::array< std::array< std::unique_ptr< TH1 >, resolTypes >, samplingsm_hConst
const TArrayD * m_etaBins

Detailed Description

get resolution for electron and photons (converted / unconverted) vs E,eta

Different parameterizations (gaussian core, sigma eff 90% and 80% from crystal ball fits) are available.

This in MC based without pileup.

Definition at line 27 of file eg_resolution.h.

Constructor & Destructor Documentation

◆ eg_resolution()

eg_resolution::eg_resolution ( const std::string & configuration)

constructor (initialization done there reading root files with resolution fit parameters

Definition at line 29 of file eg_resolution.cxx.

30{
31
32 static const std::map<std::string, std::array<std::string, 4>> fileMap = {
33 {"run1",
34 {"ElectronPhotonFourMomentumCorrection/v5/"
35 "resolutionFit_electron_run1.root",
36 "ElectronPhotonFourMomentumCorrection/v5/"
37 "resolutionFit_recoUnconv_run1.root",
38 "ElectronPhotonFourMomentumCorrection/v5/"
39 "resolutionFit_recoConv_run1.root",
40 "ElectronPhotonFourMomentumCorrection/v5/"
41 "resolutionFit_trueUnconv_run1.root"}},
42 {"run2_pre",
43 {"ElectronPhotonFourMomentumCorrection/v5/"
44 "resolutionFit_electron_run2_pre.root",
45 "ElectronPhotonFourMomentumCorrection/v5/"
46 "resolutionFit_recoUnconv_run2_pre.root",
47 "ElectronPhotonFourMomentumCorrection/v5/"
48 "resolutionFit_recoConv_run2_pre.root",
49 "ElectronPhotonFourMomentumCorrection/v5/"
50 "resolutionFit_trueUnconv_run2_pre.root"}},
51 {"run2_R21_v1",
52 {"ElectronPhotonFourMomentumCorrection/v20/"
53 "resolutionFit_electron_run2_release21_es2017_R21_v1.root",
54 "ElectronPhotonFourMomentumCorrection/v20/"
55 "resolutionFit_recoUnconv_run2_release21_es2017_R21_v1.root",
56 "ElectronPhotonFourMomentumCorrection/v20/"
57 "resolutionFit_recoConv_run2_release21_es2017_R21_v1.root",
58 "ElectronPhotonFourMomentumCorrection/v20/"
59 "resolutionFit_trueUnconvertedPhoton_run2_release21_es2017_R21_v1.root"}}};
60 std::unique_ptr<TFile> file0;
61 std::unique_ptr<TFile> file1;
62 std::unique_ptr<TFile> file2;
63 std::unique_ptr<TFile> file3;
64
65 auto config = fileMap.find(configuration);
66
67 if (config != fileMap.end()) {
68 const std::array<std::string, 4>& files = config->second;
69 file0 =
70 std::make_unique<TFile>(PathResolverFindCalibFile(files[0]).c_str());
71 file1 =
72 std::make_unique<TFile>(PathResolverFindCalibFile(files[1]).c_str());
73 file2 =
74 std::make_unique<TFile>(PathResolverFindCalibFile(files[2]).c_str());
75 file3 =
76 std::make_unique<TFile>(PathResolverFindCalibFile(files[3]).c_str());
77 }
78
79 if (!file0 or !file1 or !file2 or !file3) {
80 throw std::runtime_error("cannot find input file for resolutions");
81 }
82
83 // samping term
84 m_hSampling[0][0] = get_object<TH1>(*file0, "hsamplingG");
85 m_hSampling[0][1] = get_object<TH1>(*file0, "hsampling80");
86 m_hSampling[0][2] = get_object<TH1>(*file0, "hsampling90");
87 m_hSampling[1][0] = get_object<TH1>(*file1, "hsamplingG");
88 m_hSampling[1][1] = get_object<TH1>(*file1, "hsampling80");
89 m_hSampling[1][2] = get_object<TH1>(*file1, "hsampling90");
90 m_hSampling[2][0] = get_object<TH1>(*file2, "hsamplingG");
91 m_hSampling[2][1] = get_object<TH1>(*file2, "hsampling80");
92 m_hSampling[2][2] = get_object<TH1>(*file2, "hsampling90");
93 m_hSampling[3][0] = get_object<TH1>(*file3, "hsamplingG");
94 m_hSampling[3][1] = get_object<TH1>(*file3, "hsampling80");
95 m_hSampling[3][2] = get_object<TH1>(*file3, "hsampling90");
96
97 // noise term
98 m_hNoise[0][0] = get_object<TH1>(*file0, "hnoiseG");
99 m_hNoise[0][1] = get_object<TH1>(*file0, "hnoise80");
100 m_hNoise[0][2] = get_object<TH1>(*file0, "hnoise90");
101 m_hNoise[1][0] = get_object<TH1>(*file1, "hnoiseG");
102 m_hNoise[1][1] = get_object<TH1>(*file1, "hnoise80");
103 m_hNoise[1][2] = get_object<TH1>(*file1, "hnoise90");
104 m_hNoise[2][0] = get_object<TH1>(*file2, "hnoiseG");
105 m_hNoise[2][1] = get_object<TH1>(*file2, "hnoise80");
106 m_hNoise[2][2] = get_object<TH1>(*file2, "hnoise90");
107 m_hNoise[3][0] = get_object<TH1>(*file3, "hnoiseG");
108 m_hNoise[3][1] = get_object<TH1>(*file3, "hnoise80");
109 m_hNoise[3][2] = get_object<TH1>(*file3, "hnoise90");
110
111 // constant term
112 m_hConst[0][0] = get_object<TH1>(*file0, "hconstG");
113 m_hConst[0][1] = get_object<TH1>(*file0, "hconst80");
114 m_hConst[0][2] = get_object<TH1>(*file0, "hconst90");
115 m_hConst[1][0] = get_object<TH1>(*file1, "hconstG");
116 m_hConst[1][1] = get_object<TH1>(*file1, "hconst80");
117 m_hConst[1][2] = get_object<TH1>(*file1, "hconst90");
118 m_hConst[2][0] = get_object<TH1>(*file2, "hconstG");
119 m_hConst[2][1] = get_object<TH1>(*file2, "hconst80");
120 m_hConst[2][2] = get_object<TH1>(*file2, "hconst90");
121 m_hConst[3][0] = get_object<TH1>(*file3, "hconstG");
122 m_hConst[3][1] = get_object<TH1>(*file3, "hconst80");
123 m_hConst[3][2] = get_object<TH1>(*file3, "hconst90");
124
125 TAxis* aa = m_hSampling[0][0]->GetXaxis();
126 m_etaBins = aa->GetXbins();
127}
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
std::array< std::array< std::unique_ptr< TH1 >, resolTypes >, samplings > m_hNoise
std::array< std::array< std::unique_ptr< TH1 >, resolTypes >, samplings > m_hConst
std::array< std::array< std::unique_ptr< TH1 >, resolTypes >, samplings > m_hSampling
const TArrayD * m_etaBins
std::vector< std::string > files
file names and file pointers
Definition hcg.cxx:50

◆ ~eg_resolution()

eg_resolution::~eg_resolution ( )
default

Member Function Documentation

◆ getResolution() [1/2]

double eg_resolution::getResolution ( const xAOD::Egamma & particle,
int resol_type = 2 ) const

get relative resolution (sigmaE/E) for egamma particles

resolution type : 0=gaussian core, 1=sigma eff 80%, 2=sigma eff 90% (default)

Definition at line 180 of file eg_resolution.cxx.

182{
183 int particle_type = -1;
184 if (particle.type() == xAOD::Type::Electron) {
185 particle_type = 0;
186 } else if (particle.type() == xAOD::Type::Photon) {
187 const xAOD::Photon* ph = static_cast<const xAOD::Photon*>(&particle);
188 const xAOD::Vertex* phVertex = ph->vertex();
189 if (phVertex) {
190 const Amg::Vector3D& pos = phVertex->position();
191 const double Rconv = static_cast<float>(hypot(pos.x(), pos.y()));
192 if (Rconv > 0 and Rconv < 800) {
193 particle_type = 2;
194 } else {
195 particle_type = 1;
196 }
197 } else {
198 particle_type = 1;
199 }
200 }
201 assert(particle_type != 1);
202
203 const double eta = particle.caloCluster()->eta();
204 const double energy = particle.e();
205 return getResolution(particle_type, energy, eta, resolution_type);
206}
Scalar eta() const
pseudorapidity method
double getResolution(int particle_type, double energy, double eta, int resol_type=2) const
get relative resolution (sigmaE/E) as a function of E (in Mev) and eta
const xAOD::Vertex * vertex(size_t index=0) const
Pointer to the xAOD::Vertex/es that match the photon candidate.
Definition Photon_v1.cxx:61
const Amg::Vector3D & position() const
Returns the 3-pos.
Eigen::Matrix< double, 3, 1 > Vector3D
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
@ Photon
The object is a photon.
Definition ObjectType.h:47
@ Electron
The object is an electron.
Definition ObjectType.h:46
Vertex_v1 Vertex
Define the latest version of the vertex class.
Photon_v1 Photon
Definition of the current "egamma version".

◆ getResolution() [2/2]

double eg_resolution::getResolution ( int particle_type,
double energy,
double eta,
int resol_type = 2 ) const

get relative resolution (sigmaE/E) as a function of E (in Mev) and eta

particle type : 0=electron, 1=reco unconverted photon, 2=reco converted photon, 3=true unconverted photon

resolution type : 0=gaussian core, 1=sigma eff 80%, 2=sigma eff 90% (default)

Definition at line 136 of file eg_resolution.cxx.

140{
141
142 if (particle_type < 0 || particle_type > 3) {
143 throw std::runtime_error("particle type must be 1, 2 or 3");
144 }
145
146 if (resolution_type < 0 || resolution_type > 2) {
147 throw std::runtime_error("resolution type must be 0, 1, 2");
148 }
149
150 const float aeta = std::abs(eta);
151 //For eta outside parameterisation range, uses last bin
152 int ibinEta = m_etaBins->GetSize() - 1;
153 for (int i = 1; i < m_etaBins->GetSize(); ++i) {
154 if (aeta < m_etaBins->GetAt(i)) {
155 ibinEta = i;
156 break;
157 }
158 }
159
160 // This can never happen ! To remove
161 if (ibinEta < 0 || ibinEta >= m_etaBins->GetSize()) {
162 throw std::runtime_error("eta outside range");
163 }
164
165 const double energyGeV = energy * 1E-3;
166 const double rsampling =
167 m_hSampling[particle_type][resolution_type]->GetBinContent(ibinEta);
168 const double rnoise =
169 m_hNoise[particle_type][resolution_type]->GetBinContent(ibinEta);
170 const double rconst =
171 m_hConst[particle_type][resolution_type]->GetBinContent(ibinEta);
172
173 const double sigma2 = rsampling * rsampling / energyGeV +
174 rnoise * rnoise / energyGeV / energyGeV +
175 rconst * rconst;
176 return sqrt(sigma2);
177}

Member Data Documentation

◆ m_etaBins

const TArrayD* eg_resolution::m_etaBins
private

Definition at line 60 of file eg_resolution.h.

◆ m_hConst

std::array<std::array<std::unique_ptr<TH1>, resolTypes>, samplings> eg_resolution::m_hConst
private

Definition at line 59 of file eg_resolution.h.

◆ m_hNoise

std::array<std::array<std::unique_ptr<TH1>, resolTypes>, samplings> eg_resolution::m_hNoise
private

Definition at line 58 of file eg_resolution.h.

◆ m_hSampling

std::array<std::array<std::unique_ptr<TH1>, resolTypes>, samplings> eg_resolution::m_hSampling
private

Definition at line 57 of file eg_resolution.h.

◆ resolTypes

size_t eg_resolution::resolTypes = 3
staticconstexpr

Definition at line 32 of file eg_resolution.h.

◆ samplings

size_t eg_resolution::samplings = 4
staticconstexpr

Definition at line 31 of file eg_resolution.h.


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