ATLAS Offline Software
Loading...
Searching...
No Matches
TFCSBinnedShowerBase.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include "HepPDT/ParticleData.hh" //for HepPDT::ParticleID
11
12#include "TMath.h"
13
14#include <cmath>
15#include <cstdlib>
16#include <limits>
17
18//=============================================
19//======= TFCSBinnedShowerBase =========
20//=============================================
21
22TFCSBinnedShowerBase::TFCSBinnedShowerBase(const char *name, const char *title)
25}
26
28
30 TFCSSimulationState &simulstate, const TFCSTruthState *truth,
31 const TFCSExtrapolationState *extrapol) const {
32
33 // select a random event from the library
34 float eta_center, phi_center;
35 long unsigned int reference_layer_index = CaloCell_ID_FCS::CaloSample_FCS::EMB2;
36 eta_center =
37 extrapol->eta(reference_layer_index, TFCSExtrapolationState::SUBPOS_MID);
38 if (eta_center > 1.4) { // Endcap becomes more relevant
39 reference_layer_index = CaloCell_ID_FCS::CaloSample_FCS::EME2;
40 }
41 // TODO: What about the endcap?
42 eta_center =
43 extrapol->eta(reference_layer_index, TFCSExtrapolationState::SUBPOS_MID);
44 phi_center =
45 extrapol->phi(reference_layer_index, TFCSExtrapolationState::SUBPOS_MID);
46
47 // Fill the total energy and layer energies into simulstate
48 float Einit;
49 const float Ekin = truth->Ekin();
50
51 if (OnlyScaleEnergy())
52 Einit = simulstate.E();
53 else
54 Einit = Ekin;
55
56 // Reset the total energy
57 simulstate.set_E(0);
58
59 get_event(simulstate, eta_center, phi_center, Einit, reference_layer_index);
60
61 for (long unsigned int layer_index = 0;
62 layer_index < CaloCell_ID_FCS::MaxSample; ++layer_index) {
63
64 float layer_energy = get_layer_energy(simulstate, layer_index);
65
66 // Reset and set the layer energy
67 simulstate.set_E(layer_index, 0);
68 simulstate.add_E(layer_index, layer_energy);
69 }
70
71 if (simulstate.E() > std::numeric_limits<double>::epsilon()) {
72 for (int ilayer = 0; ilayer < CaloCell_ID_FCS::MaxSample; ++ilayer) {
73 simulstate.set_Efrac(ilayer, simulstate.E(ilayer) / simulstate.E());
74 }
75 }
76 return FCSSuccess;
77}
78
80 Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState *truth,
81 const TFCSExtrapolationState */*extrapol*/) {
82
83 const int pdgId = truth->pdgid();
84 const float charge = HepPDT::ParticleID(pdgId).charge();
85 long unsigned int layer_index = calosample();
86
87 const double center_eta = hit.center_eta();
88 const double center_phi = hit.center_phi();
89 const double center_r = hit.center_r();
90 const double center_z = hit.center_z();
91
92 ATH_MSG_VERBOSE(" Layer " << layer_index << " Extrap eta " << center_eta
93 << " phi " << center_phi << " R " << center_r);
94
95 //next MR: change to std::abs, std::sqrt functions
96 const float dist000 = TMath::Sqrt(center_r * center_r + center_z * center_z);
97 const float eta_jakobi = TMath::Abs(2.0 * TMath::Exp(-center_eta) /
98 (1.0 + TMath::Exp(-2 * center_eta)));
99
100 long unsigned int hit_index = hit.idx();
101
102 // Get necessary the hit information
103 float r, alpha, E;
104 std::tie(r, alpha, E) =
105 get_hit_position_and_energy(simulstate, layer_index, hit_index);
106
107 hit.reset();
108 hit.E() = E;
109
110 if (layer_index <= CaloCell_ID_FCS::CaloSample_FCS::FCAL0) {
111 float delta_eta_mm = r * std::cos(alpha);
112 float delta_phi_mm = r * std::sin(alpha);
113
114 // Particles with negative eta are expected to have the same shape
115 // as those with positive eta after transformation: delta_eta -->
116 // -delta_eta
117 if (center_eta < 0.) {
118 delta_eta_mm = -delta_eta_mm;
119 }
120
121 // We derive the shower shapes for electrons and positively charged
122 // hadrons. Particle with the opposite charge are expected to have the
123 // same shower shape after the transformation: delta_phi -->
124 // -delta_phi
125 if ((charge < 0. && pdgId != 11) || pdgId == -11)
126 delta_phi_mm = -delta_phi_mm;
127
128 const float delta_eta = delta_eta_mm / eta_jakobi / dist000;
129 const float delta_phi = delta_phi_mm / center_r;
130
131 hit.eta() = center_eta + delta_eta;
132 hit.phi() = TVector2::Phi_mpi_pi(center_phi + delta_phi);
133
134 ATH_MSG_VERBOSE(" Hit eta " << hit.eta() << " phi " << hit.phi()
135 << " layer " << layer_index);
136
137 } else { // FCAL is in (x,y,z)
138 const float hit_r = r * std::cos(alpha) + center_r;
139 float delta_phi = r * std::sin(alpha) / center_r;
140 // We derive the shower shapes for electrons and positively charged
141 // hadrons. Particle with the opposite charge are expected to have the
142 // same shower shape after the transformation: delta_phi -->
143 // -delta_phi
144 if ((charge < 0. && pdgId != 11) || pdgId == -11)
145 delta_phi = -delta_phi;
146 const float hit_phi = TVector2::Phi_mpi_pi(center_phi + delta_phi);
147 hit.x() = hit_r * std::cos(hit_phi);
148 hit.y() = hit_r * std::sin(hit_phi);
149 hit.z() = center_z;
150 ATH_MSG_VERBOSE(" Hit x " << hit.x() << " y " << hit.y() << " layer "
151 << layer_index);
152 }
153
154 return FCSSuccess;
155}
#define ATH_MSG_VERBOSE(x)
double charge(const T &p)
Definition AtlasPID.h:997
FCSReturnCode
Base class for all FastCaloSim parametrizations Functionality in derivde classes is provided through ...
TFCSBinnedShowerBase(const char *name=nullptr, const char *title=nullptr)
virtual float get_layer_energy(TFCSSimulationState &simulstate, long unsigned int layer_index) const =0
virtual FCSReturnCode simulate_hit(Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) override
simulated one hit position with some energy.
virtual std::tuple< float, float, float > get_hit_position_and_energy(TFCSSimulationState &simulstate, long unsigned int layer_index, long unsigned int hit_index) const =0
virtual void get_event(TFCSSimulationState &simulstate, float eta_center, float phi_center, float e_init, long unsigned int reference_layer_index) const =0
do not persistify
virtual FCSReturnCode simulate(TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) const override
Method in all derived classes to do some simulation.
double phi(int layer, int subpos) const
double eta(int layer, int subpos) const
TFCSLateralShapeParametrizationHitBase(const char *name=nullptr, const char *title=nullptr)
void set_E(int sample, double Esample)
void add_E(int sample, double Esample)
void set_Efrac(int sample, double Efracsample)
int pdgid() const
double Ekin() const
int r
Definition globals.cxx:22