ATLAS Offline Software
TFCS2DFunctionLateralShapeParametrization.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 #include "CLHEP/Random/RandFlat.h"
6 #include "CLHEP/Random/RandPoisson.h"
7 
12 
13 #include "TFile.h"
14 #include "TMath.h"
15 #include "TH2.h"
16 
17 #include "HepPDT/ParticleData.hh"
18 #include "HepPDT/ParticleDataTable.hh"
19 
20 //=============================================
21 //======= TFCS2DFunctionLateralShapeParametrization =========
22 //=============================================
23 
26  const char *title)
27  : TFCSLateralShapeParametrizationHitBase(name, title), m_function(nullptr),
28  m_nhits(0) {
30 }
31 
34  if (m_function)
35  delete m_function;
36  m_function = nullptr;
37 }
38 
40  TFCSSimulationState & /*simulstate*/, const TFCSTruthState * /*truth*/,
41  const TFCSExtrapolationState * /*extrapol*/) const {
42  return 1.0 / m_nhits;
43 }
44 
46  TFCSSimulationState &simulstate, const TFCSTruthState * /*truth*/,
47  const TFCSExtrapolationState * /*extrapol*/) const {
48  if (!simulstate.randomEngine()) {
49  return -1;
50  }
51 
52  return CLHEP::RandPoisson::shoot(simulstate.randomEngine(), m_nhits);
53 }
54 
56  float nhits) {
57  m_nhits = nhits;
58 }
59 
61  Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState *truth,
62  const TFCSExtrapolationState * /*extrapol*/) {
63  if (!simulstate.randomEngine()) {
64  return FCSFatal;
65  }
66  if (!m_function) {
67  return FCSFatal;
68  }
69 
70  const int pdgId = truth->pdgid();
71  const double charge = HepPDT::ParticleID(pdgId).charge();
72 
73  const int cs = calosample();
74  const double center_eta = hit.center_eta();
75  const double center_phi = hit.center_phi();
76  const double center_r = hit.center_r();
77  const double center_z = hit.center_z();
78 
79  if (TMath::IsNaN(center_r) or TMath::IsNaN(center_z) or
80  TMath::IsNaN(center_eta) or
81  TMath::IsNaN(center_phi)) { // Check if extrapolation fails
82  return FCSFatal;
83  }
84 
85  float alpha, r, rnd1, rnd2;
86  rnd1 = CLHEP::RandFlat::shoot(simulstate.randomEngine());
87  rnd2 = CLHEP::RandFlat::shoot(simulstate.randomEngine());
88  if (is_phi_symmetric()) {
89  if (rnd2 >= 0.5) { // Fill negative phi half of shape
90  rnd2 -= 0.5;
91  rnd2 *= 2;
92  m_function->rnd_to_fct(alpha, r, rnd1, rnd2);
93  alpha = -alpha;
94  } else { // Fill positive phi half of shape
95  rnd2 *= 2;
96  m_function->rnd_to_fct(alpha, r, rnd1, rnd2);
97  }
98  } else {
99  m_function->rnd_to_fct(alpha, r, rnd1, rnd2);
100  }
101  if (TMath::IsNaN(alpha) || TMath::IsNaN(r)) {
102  ATH_MSG_ERROR(" 2D function, #hits=" << m_nhits << " alpha=" << alpha
103  << " r=" << r << " rnd1=" << rnd1
104  << " rnd2=" << rnd2);
105  alpha = 0;
106  r = 0.001;
107 
108  ATH_MSG_ERROR(" This error could probably be retried");
109  return FCSFatal;
110  }
111 
112  float delta_eta_mm = r * cos(alpha);
113  float delta_phi_mm = r * sin(alpha);
114 
115  // Particles with negative eta are expected to have the same shape as those
116  // with positive eta after transformation: delta_eta --> -delta_eta
117  if (center_eta < 0.)
118  delta_eta_mm = -delta_eta_mm;
119  // We derive the shower shapes for electrons and positively charged hadrons.
120  // Particle with the opposite charge are expected to have the same shower shape
121  // after the transformation: delta_phi --> -delta_phi
122  if ((charge < 0. && pdgId!=11) || pdgId==-11)
123  delta_phi_mm = -delta_phi_mm;
124 
125  const float dist000 = TMath::Sqrt(center_r * center_r + center_z * center_z);
126  const float eta_jakobi = TMath::Abs(2.0 * TMath::Exp(-center_eta) /
127  (1.0 + TMath::Exp(-2 * center_eta)));
128 
129  const float delta_eta = delta_eta_mm / eta_jakobi / dist000;
130  const float delta_phi = delta_phi_mm / center_r;
131 
132  hit.setEtaPhiZE(center_eta + delta_eta, center_phi + delta_phi, center_z,
133  hit.E());
134 
135  ATH_MSG_DEBUG("HIT: E=" << hit.E() << " cs=" << cs << " eta=" << hit.eta()
136  << " phi=" << hit.phi() << " z=" << hit.z()
137  << " r=" << r << " alpha=" << alpha);
138 
139  return FCSSuccess;
140 }
141 
143  float nhits) {
144  if (!func)
145  return false;
146  if (m_function)
147  delete m_function;
148  m_function = func;
149 
150  if (nhits > 0)
151  set_number_of_hits(nhits);
152 
153  return true;
154 }
155 
157  TString opt(option);
158  bool shortprint = opt.Index("short") >= 0;
159  bool longprint = msgLvl(MSG::DEBUG) || (msgLvl(MSG::INFO) && !shortprint);
160  TString optprint = opt;
161  optprint.ReplaceAll("short", "");
163 
164  if (longprint) {
165  if (is_phi_symmetric()) {
166  ATH_MSG_INFO(optprint << " 2D function, #hits=" << m_nhits
167  << " (phi symmetric)");
168  } else {
169  ATH_MSG_INFO(optprint << " 2D function, #hits=" << m_nhits
170  << " (not phi symmetric)");
171  }
172  }
173 }
beamspotman.r
def r
Definition: beamspotman.py:676
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
TFCS2DFunctionLateralShapeParametrization::simulate_hit
virtual FCSReturnCode simulate_hit(Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) override
simulated one hit position with weight that should be put into simulstate sometime later all hit weig...
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:60
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
make_coralServer_rep.opt
opt
Definition: make_coralServer_rep.py:19
TFCS2DFunctionLateralShapeParametrization::get_sigma2_fluctuation
virtual double get_sigma2_fluctuation(TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) const override
default for this class is to simulate get_number_of_expected_hits() hits, which gives fluctuations si...
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:39
TFCSLateralShapeParametrization::Print
void Print(Option_t *option="") const override
Definition: TFCSLateralShapeParametrization.cxx:53
TFCSExtrapolationState
Definition: TFCSExtrapolationState.h:13
TFCSSimulationState::randomEngine
CLHEP::HepRandomEngine * randomEngine()
Definition: TFCSSimulationState.h:36
TFCSLateralShapeParametrizationHitBase::Hit
Definition: TFCSLateralShapeParametrizationHitBase.h:42
drawFromPickle.cos
cos
Definition: drawFromPickle.py:36
TFCSLateralShapeParametrizationHitBase::Hit::center_phi
float & center_phi()
Definition: TFCSLateralShapeParametrizationHitBase.h:101
TFCSLateralShapeParametrizationHitBase
Definition: TFCSLateralShapeParametrizationHitBase.h:13
TFCS2DFunctionLateralShapeParametrization::Print
virtual void Print(Option_t *option="") const override
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:156
TFCSLateralShapeParametrizationHitBase::Hit::E
float & E()
Definition: TFCSLateralShapeParametrizationHitBase.h:90
TFCS2DFunction::rnd_to_fct
virtual void rnd_to_fct(float &valuex, float &valuey, float rnd0, float rnd1) const =0
TFCS2DFunctionLateralShapeParametrization::TFCS2DFunctionLateralShapeParametrization
TFCS2DFunctionLateralShapeParametrization(const char *name=nullptr, const char *title=nullptr)
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:25
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
TFCS2DFunctionLateralShapeParametrization::~TFCS2DFunctionLateralShapeParametrization
~TFCS2DFunctionLateralShapeParametrization()
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:33
TFCS2DFunctionLateralShapeParametrization::get_number_of_hits
virtual int get_number_of_hits(TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) const override
default for this class is to simulate get_number_of_expected_hits() hits
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:45
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
covarianceTool.title
title
Definition: covarianceTool.py:542
TFCSLateralShapeParametrizationHitBase::Hit::center_z
float & center_z()
Definition: TFCSLateralShapeParametrizationHitBase.h:99
TFCS2DFunctionLateralShapeParametrization::reset_phi_symmetric
virtual void reset_phi_symmetric()
Definition: TFCS2DFunctionLateralShapeParametrization.h:29
jobOptions.ParticleID
ParticleID
Definition: jobOptions.decayer.py:85
FCSFatal
@ FCSFatal
Definition: TFCSParametrizationBase.h:41
TFCS2DFunctionLateralShapeParametrization::Initialize
bool Initialize(TFCS2DFunction *func, float nhits=-1)
Init from function.
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:142
FCSSuccess
@ FCSSuccess
Definition: TFCSParametrizationBase.h:41
TFCS2DFunctionLateralShapeParametrization::set_number_of_hits
void set_number_of_hits(float nhits)
set the desired number of hits
Definition: TFCS2DFunctionLateralShapeParametrization.cxx:55
TFCS2DFunctionLateralShapeParametrization.h
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
charge
double charge(const T &p)
Definition: AtlasPID.h:756
TFCSTruthState::pdgid
int pdgid() const
Definition: TFCSTruthState.h:25
TFCSLateralShapeParametrizationHitBase::Hit::z
float & z()
Definition: TFCSLateralShapeParametrizationHitBase.h:91
TFCSLateralShapeParametrization::calosample
int calosample() const
Definition: TFCSLateralShapeParametrization.h:34
TFCSExtrapolationState.h
eFEXNTuple.delta_phi
def delta_phi(phi1, phi2)
Definition: eFEXNTuple.py:15
DEBUG
#define DEBUG
Definition: page_access.h:11
TFCSLateralShapeParametrizationHitBase::Hit::center_eta
float & center_eta()
Definition: TFCSLateralShapeParametrizationHitBase.h:100
TFCS2DFunction
Definition: TFCS2DFunction.h:15
TFCSLateralShapeParametrizationHitBase::Hit::center_r
float & center_r()
Definition: TFCSLateralShapeParametrizationHitBase.h:98
TFCSLateralShapeParametrizationHitBase::Hit::eta
float & eta()
Definition: TFCSLateralShapeParametrizationHitBase.h:86
TFCSSimulationState.h
FastCaloSim_CaloCell_ID.h
drawFromPickle.sin
sin
Definition: drawFromPickle.py:36
TFCSLateralShapeParametrizationHitBase::Hit::setEtaPhiZE
void setEtaPhiZE(float eta, float phi, float z, float E)
Definition: TFCSLateralShapeParametrizationHitBase.h:56
TFCS2DFunctionLateralShapeParametrization::is_phi_symmetric
bool is_phi_symmetric() const
Definition: TFCS2DFunctionLateralShapeParametrization.h:27
TFCSTruthState
Definition: TFCSTruthState.h:13
TFCS2DFunctionLateralShapeParametrization::m_function
TFCS2DFunction * m_function
Histogram to be used for the shape simulation.
Definition: TFCS2DFunctionLateralShapeParametrization.h:67
TFCSSimulationState
Definition: TFCSSimulationState.h:32
ISF_FCS::MLogging::msgLvl
bool msgLvl(const MSG::Level lvl) const
Check whether the logging system is active at the provided verbosity level.
Definition: MLogging.h:222
TFCS2DFunctionLateralShapeParametrization::m_nhits
float m_nhits
Definition: TFCS2DFunctionLateralShapeParametrization.h:68