ATLAS Offline Software
TFCSLateralShapeTuning.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
7 
8 //=============================================
9 //======= TFCSLateralShapeTuning ==============
10 //=============================================
11 
13  const char *title)
15 
17  // clear parameter interpolation map and free memory
18  for (auto &p : m_parameterInterpol) {
19  delete p.second;
20  }
21  m_parameterInterpol.clear();
22 }
23 
25  const std::string &pathToModelParameters, int intMinEta, int intMaxEta) {
26  // get current calo layer
28 
29  // load file containing model parameters
30  std::unique_ptr<TFile> modelParametersFile = std::unique_ptr<TFile>(
31  TFile::Open(pathToModelParameters.c_str(), "READ"));
32  modelParametersFile->cd();
33 
34  // set parameter model names depending on layer
35  std::vector<std::string> parameterModelNames;
36  if (layer == 1 || layer == 5)
37  parameterModelNames = {"a0", "a1", "a2", "a3"};
38  else if (layer == 2 || layer == 6)
39  parameterModelNames = {"eta_s", "phi_s"};
40 
41  for (auto &parameterName : parameterModelNames) {
42 
43  // interpolate the model parameters
44  TFCSEnergyInterpolationPiecewiseLinear *linModelInterpol =
46 
47  TGraph *modelParameterGraph = (TGraph *)gDirectory->Get(
48  Form("photons/layer%d/m%d_m%d_%d_%d/%s", layer, intMaxEta, intMinEta,
49  intMinEta, intMaxEta, parameterName.c_str()));
50  if (modelParameterGraph) {
51  // initialize the model interpolation and save in interpolation map
52  linModelInterpol->InitFromArrayInEkin(modelParameterGraph->GetN(),
53  modelParameterGraph->GetX(),
54  modelParameterGraph->GetY());
55  m_parameterInterpol.insert(
56  std::make_pair(parameterName, linModelInterpol));
57  } else {
58  ATH_MSG_DEBUG("[TFCSLateralShapeTuning] Could not find model parameter "
59  "graph for layer="
60  << layer << " minEta=" << intMinEta
61  << " maxEta=" << intMaxEta);
62  return FCSSuccess;
63  }
64  }
65  modelParametersFile->Close();
66 
67  return FCSSuccess;
68 }
69 
72  ATH_MSG_DEBUG("[TFCSLateralShapeTuning] Initializing data tuning model from "
73  "interpolation map.");
75  return FCSSuccess;
76 }
77 
80  const TFCSTruthState *truth,
81  const TFCSExtrapolationState *) {
82 
83  // do not do anything if the parameter interpolation map is empty
84  // this means we are in an pseudorapidity region, where no tuning to data is
85  // available
86  if (m_parameterInterpol.empty())
87  return FCSSuccess;
88 
89  // set maximum scaling factor
90  float maxScaling = 500;
91 
92  // retrieve particle data
93  const int pdgId = truth->pdgid();
94  const double charge = HepPDT::ParticleID(pdgId).charge();
95 
96  // retreive hit information
97  const double centerEta = hit.center_eta();
98  const double centerPhi = hit.center_phi();
99  const double centerZ = hit.center_z();
100  // retrieve truth kinetic energy for interpolation
101  const float Ekin = truth->Ekin();
102  // retrieve calo sample
104 
105  ATH_MSG_DEBUG("[TFCSLateralShapeTuning] Initializing with pdgId="
106  << pdgId << ", charge=" << charge << ", Ekin=" << Ekin
107  << ", caloSample=" << layer);
108 
109  ATH_MSG_DEBUG("[TFCSLateralShapeTuning] Got hit position: "
110  << " hit.eta=" << hit.eta() << ", hit.phi=" << hit.phi()
111  << ", hit.r=" << hit.r());
112 
113  // compute deltaEta and deltaPhi
114  const double deltaEta = hit.eta() - centerEta;
115  const double deltaPhi = hit.phi() - centerPhi;
116 
117  if (layer == 2 || layer == 6) {
118 
119  double etaScaleFactor = m_parameterInterpol["eta_s"]->evaluate(Ekin);
120  double phiScaleFactor = m_parameterInterpol["phi_s"]->evaluate(Ekin);
121 
122  // add a maximum scaling threshold to prevent unreasonable scalings
123  etaScaleFactor =
124  std::abs(etaScaleFactor) < maxScaling ? etaScaleFactor : maxScaling;
125  phiScaleFactor =
126  std::abs(phiScaleFactor) < maxScaling ? phiScaleFactor : maxScaling;
127 
128  ATH_MSG_DEBUG("[TFCSLateralShapeTuning] Applying 2D eta_s - eta_phi "
129  "scaling model with eta_s="
130  << etaScaleFactor << " and phi_s=" << phiScaleFactor);
131 
132  double deltaEtaCorr = etaScaleFactor * deltaEta;
133  double deltaPhiCorr = phiScaleFactor * deltaPhi;
134 
135  hit.setEtaPhiZE(centerEta + deltaEtaCorr, centerPhi + deltaPhiCorr, centerZ,
136  hit.E());
137 
138  }
139 
140  else if (layer == 1 || layer == 5) {
141 
142  double etaScaleFactor = getSeriesScalingFactor(
143  m_parameterInterpol["a0"]->evaluate(Ekin),
144  m_parameterInterpol["a1"]->evaluate(Ekin),
145  m_parameterInterpol["a2"]->evaluate(Ekin),
146  m_parameterInterpol["a3"]->evaluate(Ekin), std::abs(deltaEta));
147 
148  // add a maximum scaling threshold to prevent unreasonable scalings
149  etaScaleFactor =
150  std::abs(etaScaleFactor) < maxScaling ? etaScaleFactor : maxScaling;
151 
152  ATH_MSG_DEBUG("[TFCSLateralShapeTuning] Applying eta_s series expansion "
153  "model with eta_s="
154  << etaScaleFactor);
155 
156  double deltaEtaCorr = etaScaleFactor * deltaEta;
157 
158  hit.setEtaPhiZE(centerEta + deltaEtaCorr, centerPhi + deltaPhi, centerZ,
159  hit.E());
160  }
161 
162  return FCSSuccess;
163 }
164 
166  double a0, double a1, double a2, double a3, double distToShowerCenter) {
167 
168  double meanDistToShowerCentre = 0.0039;
169  double scaleFactor =
170  1 + a0 + a1 * distToShowerCenter / meanDistToShowerCentre +
171  a2 * std::pow(distToShowerCenter / meanDistToShowerCentre, 2) +
172  a3 * std::pow(distToShowerCenter / meanDistToShowerCentre, 3);
173 
174  return scaleFactor;
175 }
TFCSLateralShapeTuning::m_parameterInterpol
interpolationMap m_parameterInterpol
Definition: TFCSLateralShapeTuning.h:48
FCSReturnCode
FCSReturnCode
Base class for all FastCaloSim parametrizations Functionality in derivde classes is provided through ...
Definition: TFCSParametrizationBase.h:41
TFCSLateralShapeParametrizationHitBase::Hit::phi
float & phi()
Definition: TFCSLateralShapeParametrizationHitBase.h:87
TFCSEnergyInterpolationPiecewiseLinear::InitFromArrayInEkin
void InitFromArrayInEkin(Int_t np, const Double_t Ekin[], const Double_t response[])
Definition: TFCSEnergyInterpolationPiecewiseLinear.cxx:43
TFCSLateralShapeTuning::TFCSLateralShapeTuning
TFCSLateralShapeTuning(const char *name=nullptr, const char *title=nullptr)
Constructor.
Definition: TFCSLateralShapeTuning.cxx:12
TFCSLateralShapeTuning::initFromModelFile
FCSReturnCode initFromModelFile(const std::string &pathToModelParameters, int intMinEta, int intMaxEta)
Definition: TFCSLateralShapeTuning.cxx:24
xAOD::deltaPhi
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner setEtaMap setEtaBin setIsTgcFailure setDeltaPt deltaPhi
Definition: L2StandAloneMuon_v1.cxx:160
conifer::pow
constexpr int pow(int x)
Definition: conifer.h:20
TFCSLateralShapeTuning::interpolationMap
std::map< std::string, TFCSEnergyInterpolationPiecewiseLinear * > interpolationMap
Definition: TFCSLateralShapeTuning.h:23
TFCSEnergyInterpolationPiecewiseLinear
Definition: TFCSEnergyInterpolationPiecewiseLinear.h:19
TFCSExtrapolationState
Definition: TFCSExtrapolationState.h:13
TFCSLateralShapeParametrizationHitBase::Hit
Definition: TFCSLateralShapeParametrizationHitBase.h:42
TFCSLateralShapeParametrizationHitBase::Hit::center_phi
float & center_phi()
Definition: TFCSLateralShapeParametrizationHitBase.h:101
TFCSLateralShapeTuning::initFromMap
FCSReturnCode initFromMap(const interpolationMap &)
Definition: TFCSLateralShapeTuning.cxx:71
TFCSLateralShapeTuning::~TFCSLateralShapeTuning
~TFCSLateralShapeTuning()
Definition: TFCSLateralShapeTuning.cxx:16
TFCSLateralShapeParametrizationHitBase
Definition: TFCSLateralShapeParametrizationHitBase.h:13
LArNewCalib_Delay_OFC_Cali.evaluate
evaluate
Definition: LArNewCalib_Delay_OFC_Cali.py:208
TFCSLateralShapeParametrizationHitBase::Hit::E
float & E()
Definition: TFCSLateralShapeParametrizationHitBase.h:90
TFCSLateralShapeTuning::simulate_hit
virtual FCSReturnCode simulate_hit(Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) override
simulated one hit position with some energy.
Definition: TFCSLateralShapeTuning.cxx:79
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:210
P4Helpers::deltaEta
double deltaEta(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
Definition: P4Helpers.h:66
ParseInputs.gDirectory
gDirectory
Definition: Final2012/ParseInputs.py:133
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
TFCSLateralShapeTuning.h
TFCSTruthState::Ekin
double Ekin() const
Definition: TFCSTruthState.h:26
TRT::Hit::layer
@ layer
Definition: HitInfo.h:79
covarianceTool.title
title
Definition: covarianceTool.py:542
TFCSLateralShapeParametrizationHitBase::Hit::center_z
float & center_z()
Definition: TFCSLateralShapeParametrizationHitBase.h:99
jobOptions.ParticleID
ParticleID
Definition: jobOptions.decayer.py:85
TFCSLateralShapeTuning::getSeriesScalingFactor
static double getSeriesScalingFactor(double a0, double a1, double a2, double a3, double distToShowerCenter)
Definition: TFCSLateralShapeTuning.cxx:165
a0
double a0
Definition: globals.cxx:27
FCSSuccess
@ FCSSuccess
Definition: TFCSParametrizationBase.h:41
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
charge
double charge(const T &p)
Definition: AtlasPID.h:538
TFCSTruthState::pdgid
int pdgid() const
Definition: TFCSTruthState.h:25
TFCSLateralShapeParametrizationHitBase::Hit::r
float r()
Definition: TFCSLateralShapeParametrizationHitBase.h:92
TFCSLateralShapeParametrization::calosample
int calosample() const
Definition: TFCSLateralShapeParametrization.h:34
TFCSLateralShapeParametrizationHitBase::Hit::center_eta
float & center_eta()
Definition: TFCSLateralShapeParametrizationHitBase.h:100
TFCSLateralShapeParametrizationHitBase::Hit::eta
float & eta()
Definition: TFCSLateralShapeParametrizationHitBase.h:86
TFCSSimulationState.h
TFCSLateralShapeParametrizationHitBase::Hit::setEtaPhiZE
void setEtaPhiZE(float eta, float phi, float z, float E)
Definition: TFCSLateralShapeParametrizationHitBase.h:56
TFCSTruthState
Definition: TFCSTruthState.h:13
TFCSSimulationState
Definition: TFCSSimulationState.h:32