ATLAS Offline Software
CaloSwEtamod_v2.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3 */
34 #include "CaloSwEtamod_v2.h"
37 #include <cmath>
38 
39 
40 using xAOD::CaloCluster;
42 using std::abs;
43 
44 
66  CaloCluster* cluster,
67  const CaloDetDescrElement* elt,
68  float eta,
69  float adj_eta,
70  float /*phi*/,
71  float /*adj_phi*/,
72  CaloSampling::CaloSample /*samp*/)
73  const
74 {
75  // Compute the eta offset within the cell.
76  float etamod = eta - elt->eta() + elt->deta()/2;
77  if (elt->eta_raw() < 0)
78  etamod = elt->deta() - etamod;
79 
80  // These can happen due to DD bugs.
81  if (etamod < 0)
82  etamod = 0;
83  else if (etamod > elt->deta())
84  etamod = elt->deta();
85 
86  float adj_aeta = std::abs (adj_eta);
87 
88  // ??? In most of the calorimeter, the cells in sampling 2 have an eta
89  // width of 0.025. This correction was derived with this assumption.
90  // However, it's not completely accurate: in the gap region, some
91  // cells have a different size. For example, the last sampling 2
92  // cell in the barrel (starting at 1.4) has a width of 0.075.
93  // For this pass, fix up etamod to match what was used when
94  // the correction was derived. This should be cleaned up in
95  // a subsequent version.
96  if (adj_aeta >= 1.4 && adj_aeta <= 1.475)
97  etamod = fmod (adj_aeta, 0.025);
98 
99  // Before doing the energy interpolation, make a crude total correction
100  // of the energy. This is needed since the corrections are tabulated
101  // using the true cluster energies.
102  float energy = cluster->e();
103  float rfac = interpolate (m_rfac(myctx), adj_aeta, m_rfac_degree(myctx));
104  energy /= rfac;
105 
106  // Calculate the correction.
107  float corr = energy_interpolation (energy,
108  Builder (m_correction(myctx), etamod),
109  m_energies(myctx),
110  m_energy_degree(myctx));
111 
112  // set energy, and rescale each sampling
113  setenergy (cluster, cluster->e() / corr);
114 }
115 
116 
123  float etamod)
124  : m_correction (correction),
125  m_etamod (etamod)
126 {
127 }
128 
129 
136 float CaloSwEtamod_v2::Builder::calculate (int energy_ndx, bool& good) const
137 {
138  good = true;
139  return m_correction[energy_ndx][0]
140  + m_etamod*m_correction[energy_ndx][1]
141  + m_etamod*m_etamod*m_correction[energy_ndx][2];
142 }
CaloDetDescrElement::deta
float deta() const
cell deta
Definition: Calorimeter/CaloDetDescr/CaloDetDescr/CaloDetDescrElement.h:356
CaloSwEtamod_v2::makeTheCorrection
virtual void makeTheCorrection(const Context &myctx, xAOD::CaloCluster *cluster, const CaloDetDescrElement *elt, float eta, float adj_eta, float phi, float adj_phi, CaloSampling::CaloSample samp) const override
Virtual function for the correction-specific code.
Definition: CaloSwEtamod_v2.cxx:65
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:83
CaloDetDescrElement
This class groups all DetDescr information related to a CaloCell. Provides a generic interface for al...
Definition: Calorimeter/CaloDetDescr/CaloDetDescr/CaloDetDescrElement.h:66
CaloClusterCorrection::setenergy
virtual void setenergy(xAOD::CaloCluster *cluster, float energy) const
Definition: CaloClusterCorrection.cxx:94
xAOD::CaloCluster
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
Definition: Event/xAOD/xAODCaloEvent/xAODCaloEvent/CaloCluster.h:19
CaloDetDescrManager.h
Definition of CaloDetDescrManager.
CaloDetDescrElement::eta_raw
float eta_raw() const
cell eta_raw
Definition: Calorimeter/CaloDetDescr/CaloDetDescr/CaloDetDescrElement.h:350
CaloSwCorrections.etamod
def etamod(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:206
tools.zlumi_mc_cf.correction
def correction(mu, runmode, campaign, run=None)
Definition: zlumi_mc_cf.py:4
CaloSwEtamod_v2::m_energy_degree
Constant< int > m_energy_degree
Calibration constant: degree of the polynomial interpolation in energy.
Definition: CaloSwEtamod_v2.h:161
CaloSwCorrections.rfac
def rfac(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:182
ParticleGun_FastCalo_ChargeFlip_Config.energy
energy
Definition: ParticleGun_FastCalo_ChargeFlip_Config.py:78
CaloSampling::CaloSample
CaloSample
Definition: Calorimeter/CaloGeoHelpers/CaloGeoHelpers/CaloSampling.h:22
CaloSwEtamod_v2::m_correction
Constant< CxxUtils::Array< 2 > > m_correction
Calibration constant: tabulated arrays of function parameters.
Definition: CaloSwEtamod_v2.h:142
CaloSwEtamod_v2.h
EM calorimeter modulation corrections.
CaloSwEtamod_v2::Builder::calculate
virtual float calculate(int energy_ndx, bool &good) const
Calculate the correction for tabulated energy ENERGY_NDX.
Definition: CaloSwEtamod_v2.cxx:136
CxxUtils::Array< 2 >
CaloCluster
Principal data class for CaloCell clusters.
Definition: Calorimeter/CaloEvent/CaloEvent/CaloCluster.h:79
CaloSwEtamod_v2::m_rfac_degree
Constant< int > m_rfac_degree
Definition: CaloSwEtamod_v2.h:152
CaloSwEtamod_v2::Builder::Builder
Builder(const CxxUtils::Array< 2 > &corr, float etamod)
Constructor.
Definition: CaloSwEtamod_v2.cxx:122
CaloSwEtamod_v2::Builder
Helper class for calculating the energy interpolation table.
Definition: CaloSwEtamod_v2.h:123
CaloSwEtamod_v2::m_energies
Constant< CxxUtils::Array< 1 > > m_energies
Calibration constant: table of energies at which the correction was tabulated.
Definition: CaloSwEtamod_v2.h:157
CaloSwEtamod_v2::m_rfac
Constant< CxxUtils::Array< 2 > > m_rfac
Definition: CaloSwEtamod_v2.h:147
ReadBchFromCool.good
good
Definition: ReadBchFromCool.py:433
CaloClusterCorrectionCommon::energy_interpolation
static float energy_interpolation(float energy, const TableBuilder &builder, const CaloRec::Array< 1 > &energies, int energy_degree)
Many of the corrections use the same method for doing the final interpolation in energy.
Definition: CaloClusterCorrectionCommon.cxx:575
CaloClusterCorr::interpolate
float interpolate(const CaloRec::Array< 2 > &a, float x, unsigned int degree, unsigned int ycol=1, const CaloRec::Array< 1 > &regions=CaloRec::Array< 1 >(), int n_points=-1, bool fixZero=false)
Polynomial interpolation in a table.
Definition: interpolate.cxx:75
CaloCluster::e
virtual double e() const
Retrieve energy independent of signal state.
Definition: Calorimeter/CaloEvent/CaloEvent/CaloCluster.h:753
CaloUtils::ToolConstantsContext
Context object for retrieving ToolConstant values.
Definition: ToolWithConstants.h:61
CaloDetDescrElement::eta
float eta() const
cell eta
Definition: Calorimeter/CaloDetDescr/CaloDetDescr/CaloDetDescrElement.h:344
interpolate.h
Polynomial interpolation in a table.