ATLAS Offline Software
Loading...
Searching...
No Matches
EnergyDepositionTool.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
10
11#ifndef PIXELDIGITIZATION_EnergyDepositionTool_H
12#define PIXELDIGITIZATION_EnergyDepositionTool_H
13
14#include "GaudiKernel/ServiceHandle.h"
16#include "GaudiKernel/ToolHandle.h"
17#include "Gaudi/Property.h"
18#include "GaudiKernel/Service.h"
20#include "BichselData.h"
21
24#include <vector>
25#include <utility> //for std::pair
26
27
28class PixelID;
29namespace InDetDD{
31}
32
33namespace CLHEP{
34 class HepRandomEngine;
35}
36class SiHit;
37
38
39
40//====================
41// C L A S S D E F
42//====================
44public:
45 EnergyDepositionTool(const std::string& type, const std::string& name, const IInterface* parent);
46
47 static const InterfaceID& interfaceID();
48 virtual StatusCode initialize();
49 virtual StatusCode finalize();
50
52 StatusCode initTools();
53
54 std::vector<std::pair<double, double> > bichselSim(double BetaGamma, int ParticleType, double TotalLength,
55 double InciEnergy, CLHEP::HepRandomEngine* rndmEngine) const; // output hit record in the format (hit position, energy loss)
56
57 std::vector<std::pair<double, double> > clusterHits(std::vector<std::pair<double, double> >& rawHitRecord,
58 int n_pieces) const; // cluster hits into n steps (there could be thousands of hit)
59
60 virtual StatusCode depositEnergy(const TimedHitPtr<SiHit>& phit, const InDetDD::SiDetectorElement& Module,
61 std::vector<std::pair<double, double> >& trfHitRecord,
62 std::vector<double>& initialConditions, CLHEP::HepRandomEngine* rndmEngine,
63 const EventContext &ctx) const;
64
65 // Variables
66private:
68
70 nullptr
71 };
72
73 std::vector<BichselData> m_bichselData; // vector to store Bichsel Data. Each entry is for one particle type
74
75 Gaudi::Property<int> m_numberOfSteps
76 {
77 this, "numberOfSteps", 50, "Geant4:number of steps for PixelPlanar"
78 };
79
80 Gaudi::Property<int> m_numberOfCharges
81 {
82 this, "numberOfCharges", 10, "Geant4:number of charges for PixelPlanar"
83 };
84
85 Gaudi::Property<bool> m_disableDistortions
86 {
87 this, "DisableDistortions", false, "Disable simulation of module distortions"
88 };
89
90 Gaudi::Property<bool> m_doBichsel
91 {
92 this, "doBichsel", true, "re-do charge deposition following Bichsel model"
93 };
94
95 Gaudi::Property<double> m_doBichselBetaGammaCut
96 {
97 this, "doBichselBetaGammaCut", 0.1, "minimum beta-gamma for particle to be re-simulated through Bichsel Model"
98 };
99
100 Gaudi::Property<bool> m_doDeltaRay
101 {
102 this, "doDeltaRay", false, "whether we simulate delta-ray using Bichsel model"
103 };
104
105 Gaudi::Property<double> m_DeltaRayCut
106 {
107 this, "DeltaRayCut", 80.7687, "Cut of delta ray [keV] - Value should be consistent with range cut in simulation"
108 };
109
110 Gaudi::Property<bool> m_doPU
111 {
112 this, "doPU", true, "Whether we apply Bichsel model on PU"
113 };
114
115 Gaudi::Property<int> m_nCols
116 {
117 this, "nCols", 1, "Number of collision for each sampling"
118 };
119
120 Gaudi::Property<int> m_LoopLimit
121 {
122 this, "LoopLimit", 100000, "Limit assuming 1 collision per sampling"
123 };
124
126 {
127 this, "PixelDistortionData", "PixelDistortionData", "Output readout distortion data"
128 };
129private:
130 void simulateBow(const InDetDD::SiDetectorElement* element, double& xi, double& yi, const double zi, double& xf,
131 double& yf, const double zf) const;
132 // ! unit is eV
133};
134
135#endif //PIXELDIGITIZATION_EnergyDepositionTool_H
Hold pixel distortion data produced by PixelDistortionAlg.
ParticleType
Definition TruthClasses.h:8
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Gaudi::Property< bool > m_doDeltaRay
Gaudi::Property< bool > m_doPU
Gaudi::Property< bool > m_disableDistortions
StatusCode initTools()
Gaudi::Property< double > m_doBichselBetaGammaCut
EnergyDepositionTool(const std::string &type, const std::string &name, const IInterface *parent)
Gaudi::Property< int > m_numberOfCharges
virtual StatusCode depositEnergy(const TimedHitPtr< SiHit > &phit, const InDetDD::SiDetectorElement &Module, std::vector< std::pair< double, double > > &trfHitRecord, std::vector< double > &initialConditions, CLHEP::HepRandomEngine *rndmEngine, const EventContext &ctx) const
Gaudi::Property< bool > m_doBichsel
virtual StatusCode finalize()
void simulateBow(const InDetDD::SiDetectorElement *element, double &xi, double &yi, const double zi, double &xf, double &yf, const double zf) const
std::vector< std::pair< double, double > > clusterHits(std::vector< std::pair< double, double > > &rawHitRecord, int n_pieces) const
Gaudi::Property< double > m_DeltaRayCut
Gaudi::Property< int > m_nCols
virtual ~EnergyDepositionTool()
Gaudi::Property< int > m_LoopLimit
SG::ReadCondHandleKey< PixelDistortionData > m_distortionKey
virtual StatusCode initialize()
std::vector< BichselData > m_bichselData
static const InterfaceID & interfaceID()
Gaudi::Property< int > m_numberOfSteps
std::vector< std::pair< double, double > > bichselSim(double BetaGamma, int ParticleType, double TotalLength, double InciEnergy, CLHEP::HepRandomEngine *rndmEngine) const
Class to hold geometrical description of a silicon detector element.
This is an Identifier helper class for the Pixel subdetector.
Definition PixelID.h:67
Definition SiHit.h:19
a smart pointer to a hit that also provides access to the extended timing info of the host event.
Definition TimedHitPtr.h:18
Message Stream Member.