ATLAS Offline Software
Loading...
Searching...
No Matches
TFCSFlatLateralShapeParametrization.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2019 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//=============================================
18//======= TFCSFlatLateralShapeParametrization =========
19//=============================================
20
25
27
29 TFCSSimulationState &simulstate, const TFCSTruthState * /*truth*/,
30 const TFCSExtrapolationState * /*extrapol*/) const {
31 if (!simulstate.randomEngine()) {
32 return -1;
33 }
34
35 return CLHEP::RandPoisson::shoot(simulstate.randomEngine(), m_nhits);
36}
37
41
43
45 m_scale = _scale;
46}
47
49 Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState * /*truth*/,
50 const TFCSExtrapolationState * /*extrapol*/) {
51 if (!simulstate.randomEngine()) {
52 return FCSFatal;
53 }
54
55 const int cs = calosample();
56 const double center_eta = hit.center_eta();
57 const double center_phi = hit.center_phi();
58 const double center_r = hit.center_r();
59 const double center_z = hit.center_z();
60
61 if (TMath::IsNaN(center_r) or TMath::IsNaN(center_z) or
62 TMath::IsNaN(center_eta) or
63 TMath::IsNaN(center_phi)) { // Check if extrapolation fails
64 return FCSFatal;
65 }
66
67 float alpha, r;
68
69 alpha = 2 * TMath::Pi() * CLHEP::RandFlat::shoot(simulstate.randomEngine());
70 r = m_dR * CLHEP::RandFlat::shoot(simulstate.randomEngine());
71
72 float delta_eta = r * cos(alpha);
73 float delta_phi = r * sin(alpha);
74
75 hit.setEtaPhiZE(center_eta + delta_eta, center_phi + delta_phi, center_z,
76 hit.E() * m_scale);
77
78 ATH_MSG_DEBUG("HIT: E=" << hit.E() << " cs=" << cs << " eta=" << hit.eta()
79 << " phi=" << hit.phi() << " z=" << hit.z()
80 << " r=" << r << " alpha=" << alpha);
81
82 return FCSSuccess;
83}
84
85void TFCSFlatLateralShapeParametrization::Print(Option_t *option) const {
86 TString opt(option);
87 bool shortprint = opt.Index("short") >= 0;
88 bool longprint = msgLvl(MSG::DEBUG) || (msgLvl(MSG::INFO) && !shortprint);
89 TString optprint = opt;
90 optprint.ReplaceAll("short", "");
92
93 if (longprint) {
94 ATH_MSG_INFO(optprint << " dR=" << m_dR << " scale factor=" << m_scale
95 << ", #hits=" << m_nhits);
96 }
97}
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
FCSReturnCode
Base class for all FastCaloSim parametrizations Functionality in derivde classes is provided through ...
bool msgLvl(const MSG::Level lvl) const
Check whether the logging system is active at the provided verbosity level.
Definition MLogging.h:222
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...
virtual int get_number_of_hits(TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) const override
default for this class is to simulate poisson(integral histogram) hits
void set_dR(float _dR)
set the radius in which hits should be generated
void set_number_of_hits(float nhits)
set the integral of the histogram to the desired number of hits
TFCSFlatLateralShapeParametrization(const char *name=nullptr, const char *title=nullptr)
virtual void Print(Option_t *option="") const override
void set_scale(float _scale)
set the radius in which hits should be generated
void setEtaPhiZE(float eta, float phi, float z, float E)
TFCSLateralShapeParametrizationHitBase(const char *name=nullptr, const char *title=nullptr)
void Print(Option_t *option="") const override
CLHEP::HepRandomEngine * randomEngine()
int r
Definition globals.cxx:22