ATLAS Offline Software
Loading...
Searching...
No Matches
TFCSHistoLateralShapeParametrizationFCal.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2018 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "CLHEP/Random/RandFlat.h"
6#include "CLHEP/Random/RandPoisson.h"
7
11
12#include "TMath.h"
13
14#include "HepPDT/ParticleData.hh"
15#include "HepPDT/ParticleDataTable.hh"
16
17//=============================================
18//======= TFCSHistoLateralShapeParametrizationFCal =========
19//=============================================
20
25
28
30 Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState *truth,
31 const TFCSExtrapolationState * /*extrapol*/) {
32 if (!simulstate.randomEngine()) {
33 return FCSFatal;
34 }
35
36 const int pdgId = truth->pdgid();
37 const double charge = HepPDT::ParticleID(pdgId).charge();
38
39 const int cs = calosample();
40 // const double center_phi=0.5*( extrapol->phi(cs, CaloSubPos::SUBPOS_ENT) +
41 // extrapol->phi(cs, CaloSubPos::SUBPOS_EXT) ); const double center_r=0.5*(
42 // extrapol->r(cs, CaloSubPos::SUBPOS_ENT) + extrapol->r(cs,
43 // CaloSubPos::SUBPOS_EXT) ); const double center_z=0.5*( extrapol->z(cs,
44 // CaloSubPos::SUBPOS_ENT) + extrapol->z(cs, CaloSubPos::SUBPOS_EXT) );
45 const double center_phi = hit.center_phi();
46 const double center_r = hit.center_r();
47 const double center_z = hit.center_z();
48
49 float alpha, r, rnd1, rnd2;
50 rnd1 = CLHEP::RandFlat::shoot(simulstate.randomEngine());
51 rnd2 = CLHEP::RandFlat::shoot(simulstate.randomEngine());
52 if (is_phi_symmetric()) {
53 if (rnd2 >= 0.5) { // Fill negative phi half of shape
54 rnd2 -= 0.5;
55 rnd2 *= 2;
56 m_hist.rnd_to_fct(alpha, r, rnd1, rnd2);
57 alpha = -alpha;
58 } else { // Fill positive phi half of shape
59 rnd2 *= 2;
60 m_hist.rnd_to_fct(alpha, r, rnd1, rnd2);
61 }
62 } else {
63 m_hist.rnd_to_fct(alpha, r, rnd1, rnd2);
64 }
65 if (TMath::IsNaN(alpha) || TMath::IsNaN(r)) {
66 ATH_MSG_ERROR(" Histogram: "
67 << m_hist.get_HistoBordersx().size() - 1 << "*"
68 << m_hist.get_HistoBordersy().size() - 1
69 << " bins, #hits=" << m_nhits << " alpha=" << alpha
70 << " r=" << r << " rnd1=" << rnd1 << " rnd2=" << rnd2);
71 alpha = 0;
72 r = 0.001;
73
74 ATH_MSG_ERROR(" This error could probably be retried");
75 return FCSFatal;
76 }
77
78 const float hit_r = r * cos(alpha) + center_r;
79 float delta_phi = r * sin(alpha) / center_r;
80 // We derive the shower shapes for electrons and positively charged hadrons.
81 // Particle with the opposite charge are expected to have the same shower shape
82 // after the transformation: delta_phi --> -delta_phi
83 if ((charge < 0. && pdgId!=11) || pdgId==-11)
84 delta_phi = -delta_phi;
85 const float hit_phi = delta_phi + center_phi;
86
87 hit.setXYZE(hit_r * cos(hit_phi), hit_r * sin(hit_phi), center_z, hit.E());
88
89 ATH_MSG_DEBUG("HIT: E=" << hit.E() << " cs=" << cs << " x=" << hit.x()
90 << " y=" << hit.y() << " z=" << hit.z() << " r=" << r
91 << " alpha=" << alpha);
92
93 return FCSSuccess;
94}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
double charge(const T &p)
Definition AtlasPID.h:997
FCSReturnCode
Base class for all FastCaloSim parametrizations Functionality in derivde classes is provided through ...
TFCSHistoLateralShapeParametrizationFCal(const char *name=nullptr, const char *title=nullptr)
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...
TFCS2DFunctionHistogram m_hist
Histogram to be used for the shape simulation.
TFCSHistoLateralShapeParametrization(const char *name=nullptr, const char *title=nullptr)
CLHEP::HepRandomEngine * randomEngine()
int pdgid() const
int r
Definition globals.cxx:22