ATLAS Offline Software
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 
22  const char *name, const char *title)
24  m_scale(1) {}
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 
39  m_nhits = nhits;
40 }
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 
85 void 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 }
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
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
make_coralServer_rep.opt
opt
Definition: make_coralServer_rep.py:19
TFCSLateralShapeParametrization::Print
void Print(Option_t *option="") const override
Definition: TFCSLateralShapeParametrization.cxx:53
TFCSFlatLateralShapeParametrization::set_scale
void set_scale(float _scale)
set the radius in which hits should be generated
Definition: TFCSFlatLateralShapeParametrization.cxx:44
TFCSFlatLateralShapeParametrization::m_nhits
float m_nhits
Definition: TFCSFlatLateralShapeParametrization.h:55
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
TFCSFlatLateralShapeParametrization.h
TFCSFlatLateralShapeParametrization::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 poisson(integral histogram) hits
Definition: TFCSFlatLateralShapeParametrization.cxx:28
python.utils.AtlRunQueryMemUtil._scale
dictionary _scale
Definition: AtlRunQueryMemUtil.py:7
TFCSLateralShapeParametrizationHitBase
Definition: TFCSLateralShapeParametrizationHitBase.h:13
TFCSFlatLateralShapeParametrization::m_scale
float m_scale
Definition: TFCSFlatLateralShapeParametrization.h:56
TFCSLateralShapeParametrizationHitBase::Hit::E
float & E()
Definition: TFCSLateralShapeParametrizationHitBase.h:90
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
FCSFatal
@ FCSFatal
Definition: TFCSParametrizationBase.h:41
FCSSuccess
@ FCSSuccess
Definition: TFCSParametrizationBase.h:41
TFCSFlatLateralShapeParametrization::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: TFCSFlatLateralShapeParametrization.cxx:48
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
TFCSFlatLateralShapeParametrization::~TFCSFlatLateralShapeParametrization
virtual ~TFCSFlatLateralShapeParametrization()
Definition: TFCSFlatLateralShapeParametrization.cxx:26
TFCSLateralShapeParametrizationHitBase::Hit::z
float & z()
Definition: TFCSLateralShapeParametrizationHitBase.h:91
TFCSFlatLateralShapeParametrization::Print
virtual void Print(Option_t *option="") const override
Definition: TFCSFlatLateralShapeParametrization.cxx:85
TFCSFlatLateralShapeParametrization::set_number_of_hits
void set_number_of_hits(float nhits)
set the integral of the histogram to the desired number of hits
Definition: TFCSFlatLateralShapeParametrization.cxx:38
TFCSLateralShapeParametrization::calosample
int calosample() const
Definition: TFCSLateralShapeParametrization.h:34
TFCSExtrapolationState.h
eFEXNTuple.delta_phi
def delta_phi(phi1, phi2)
Definition: eFEXNTuple.py:15
TFCSFlatLateralShapeParametrization::m_dR
float m_dR
Simulate hits flat in radius dR.
Definition: TFCSFlatLateralShapeParametrization.h:54
DEBUG
#define DEBUG
Definition: page_access.h:11
TFCSLateralShapeParametrizationHitBase::Hit::center_eta
float & center_eta()
Definition: TFCSLateralShapeParametrizationHitBase.h:100
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
TFCSFlatLateralShapeParametrization::set_dR
void set_dR(float _dR)
set the radius in which hits should be generated
Definition: TFCSFlatLateralShapeParametrization.cxx:42
TFCSTruthState
Definition: TFCSTruthState.h:13
TFCSSimulationState
Definition: TFCSSimulationState.h:32
TFCSFlatLateralShapeParametrization::TFCSFlatLateralShapeParametrization
TFCSFlatLateralShapeParametrization(const char *name=nullptr, const char *title=nullptr)
Definition: TFCSFlatLateralShapeParametrization.cxx:21
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