Loading [MathJax]/extensions/tex2jax.js
ATLAS Offline Software
Main Page
Related Pages
Modules
Namespaces
Namespace List
Namespace Members
All
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Functions
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Variables
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Typedefs
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
z
Enumerations
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
z
Enumerator
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Classes
Class List
Class Index
Class Hierarchy
Class Members
All
:
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Functions
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Variables
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Typedefs
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Enumerations
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Enumerator
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Properties
Related Functions
:
a
b
c
d
e
f
g
h
i
j
l
m
n
o
p
r
s
t
v
w
x
z
Files
File List
File Members
All
$
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Functions
_
a
b
c
d
e
f
g
h
i
j
l
m
n
o
p
q
r
s
t
u
v
w
x
z
Variables
$
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Typedefs
a
b
c
d
e
f
g
h
i
j
l
m
n
o
p
q
r
s
t
u
v
w
x
z
Enumerations
a
b
c
d
e
f
g
h
i
l
m
n
o
p
r
s
t
v
x
z
Enumerator
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Macros
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
GitLab
LXR
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Calorimeter
CaloClusterCorrection
src
CaloSwEta1e_g3.cxx
Go to the documentation of this file.
1
/*
2
Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3
*/
4
5
/********************************************************************
6
7
NAME: CaloSwEta1e_g3.cxx
8
PACKAGE: offline/Calorimeter/CaloClusterCorrection
9
10
AUTHORS: H. Ma, S. Rajagopalan
11
CREATED: Dec. 15, 1999
12
13
PURPOSE: S-shape corrections in sampling 1 as a function of eta for endcap
14
(Tuned using 100 GeV E photons)
15
base class: CaloClusterCorrection
16
Correction tuned on G3 samples.
17
18
Atrecon Orig: emreco/qeta1e.age
19
20
Updated: Feb 11, 2003 (MW)
21
new correction for DC1 data
22
Updated: Feb 28, 2003 (MW)
23
protect against unphysical clusters
24
25
Updated: May 5, 2004 (Sven Menke)
26
base class changed from algo to tool
27
28
Updated: June, 2004 (sss)
29
Use ToolWithConstants to get correction constants.
30
********************************************************************/
31
#include "
CaloSwEta1e_g3.h
"
32
#include "
CaloClusterCorrection/interpolate.h
"
33
#include <cmath>
34
35
using
xAOD::CaloCluster
;
36
using
CaloClusterCorr::interpolate
;
37
38
// granularity of middle EM barrel layer
39
40
const
float
CaloSwEta1e_g3::s_middle_layer_granularity
= 0.025;
41
42
// make correction to one cluster
43
void
CaloSwEta1e_g3::makeCorrection
(
const
Context
& myctx,
44
CaloCluster
* cluster)
const
45
{
46
// Only for endcap
47
if
(!cluster->inEndcap())
48
return
;
49
50
const
CxxUtils::Array<3>
correction
=
m_correction
(myctx);
51
const
float
correction_coef =
m_correction_coef
(myctx);
52
const
int
correction_degree =
m_correction_degree
(myctx);
53
54
float
eta
= cluster->etaSample(
CaloSampling::EME1
);
55
if
(
eta
== -999.)
return
;
56
float
etamax = cluster->etamax(
CaloSampling::EME1
);
57
float
aeta = fabs(
eta
);
58
59
int
corrndx;
60
float
ufrac;
61
if
(aeta > 2.0)
62
{
63
ufrac =
s_middle_layer_granularity
/4.;
64
corrndx = 0;
65
}
66
else
if
(aeta > 1.8)
67
{
68
ufrac =
s_middle_layer_granularity
*(1./6);
69
corrndx = 1;
70
}
71
else
if
(aeta > 1.5)
72
{
73
ufrac =
s_middle_layer_granularity
/8.;
74
corrndx = 2;
75
}
76
else
77
{
78
// not in endcap, no correction.
79
return
;
80
}
81
82
float
u
= fmod(
eta
-etamax, ufrac/2.);
83
u
+= ufrac/2.;
84
if
(
eta
< 0.)
u
= ufrac-
u
;
85
float
deta = correction_coef *
interpolate
(
correction
[corrndx],
86
u
,
87
correction_degree);
88
if
(
eta
< 0)
89
deta = -deta;
90
91
// Make the correction:
92
cluster->setEta(
CaloSampling::EME1
,
eta
- deta);
93
}
94
95
96
const
std::string&
CaloSwEta1e_g3::toolType
()
const
97
{
98
static
const
std::string
typeName
=
"CaloSwEta_g3"
;
99
return
typeName
;
100
}
eta
Scalar eta() const
pseudorapidity method
Definition:
AmgMatrixBasePlugin.h:83
CaloSwEta1e_g3::s_middle_layer_granularity
static const float s_middle_layer_granularity
Definition:
CaloSwEta1e_g3.h:83
Trk::u
@ u
Enums for curvilinear frames.
Definition:
ParamDefs.h:77
xAOD::CaloCluster
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
Definition:
Event/xAOD/xAODCaloEvent/xAODCaloEvent/CaloCluster.h:19
tools.zlumi_mc_cf.correction
def correction(mu, runmode, campaign, run=None)
Definition:
zlumi_mc_cf.py:4
CxxUtils::Array< 3 >
CaloSwEta1e_g3::makeCorrection
virtual void makeCorrection(const Context &myctx, xAOD::CaloCluster *cluster) const override
Definition:
CaloSwEta1e_g3.cxx:43
constants.EME1
int EME1
Definition:
Calorimeter/CaloClusterCorrection/python/constants.py:55
CaloCluster
Principal data class for CaloCell clusters.
Definition:
Calorimeter/CaloEvent/CaloEvent/CaloCluster.h:79
CaloSwEta1e_g3::m_correction
Constant< CxxUtils::Array< 3 > > m_correction
Definition:
CaloSwEta1e_g3.h:78
CaloSwEta1e_g3::m_correction_coef
Constant< float > m_correction_coef
Definition:
CaloSwEta1e_g3.h:79
CaloSwEta1e_g3::toolType
virtual const std::string & toolType() const override
Definition:
CaloSwEta1e_g3.cxx:96
CaloSwEta1e_g3::m_correction_degree
Constant< int > m_correction_degree
Definition:
CaloSwEta1e_g3.h:80
CaloClusterCorr::interpolate
float interpolate(const CaloRec::Array< 2 > &a, float x, unsigned int degree, unsigned int ycol=1, const CaloRec::Array< 1 > ®ions=CaloRec::Array< 1 >(), int n_points=-1, bool fixZero=false)
Polynomial interpolation in a table.
Definition:
interpolate.cxx:75
CaloUtils::ToolConstantsContext
Context object for retrieving ToolConstant values.
Definition:
ToolWithConstants.h:60
ReadCalibFromCool.typeName
typeName
Definition:
ReadCalibFromCool.py:477
interpolate.h
Polynomial interpolation in a table.
CaloSwEta1e_g3.h
Generated on Tue Mar 25 2025 21:08:13 for ATLAS Offline Software by
1.8.18