ATLAS Offline Software
Loading...
Searching...
No Matches
SCT_DetailedSurfaceChargesGenerator.h
Go to the documentation of this file.
1// -*- C++ -*-
2
3/*
4 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
5*/
6
14
15#ifndef SCT_DIGITIZATION_SCT_DETAILEDSURFACECHARGESGENERATOR_H
16#define SCT_DIGITIZATION_SCT_DETAILEDSURFACECHARGESGENERATOR_H
17
18//Inheritance
21
22// Athena
26
27// Gaudi
28#include "GaudiKernel/ITHistSvc.h" // for ITHistSvc
29#include "GaudiKernel/ServiceHandle.h"
30#include "GaudiKernel/ToolHandle.h"
31
32// CLHEP
33#include "CLHEP/Units/PhysicalConstants.h"
34
35// C++ Standard Library
36#include <string>
37
38// Charges and hits
39class SiHit;
40
41// -- do histos
42class TProfile;
43class TProfile2D;
44
45namespace InDetDD {
47}
48
49namespace CLHEP {
50 class HepRandomEngine;
51}
52
58
59class SCT_DetailedSurfaceChargesGenerator : public extends<AthAlgTool, ISurfaceChargesGenerator> {
60public:
61
63 SCT_DetailedSurfaceChargesGenerator(const std::string& type, const std::string& name, const IInterface* parent);
64
67
69 virtual StatusCode initialize() override;
70
72 virtual StatusCode finalize() override;
73
74private:
85
86 // Non-const methods are used in initialization
87 virtual void setFixedTime(float fixedTime) override { m_tfix.setValue(fixedTime); }
88
90 virtual void process(const InDetDD::SiDetectorElement* element, const TimedHitPtr<SiHit>& phit, ISiSurfaceChargesInserter& inserter, CLHEP::HepRandomEngine * rndmEngine, const EventContext& ctx) override;
91 void processSiHit(const InDetDD::SiDetectorElement* element, const SiHit& phit, ISiSurfaceChargesInserter& inserter, const float eventTime, const unsigned short eventID, CLHEP::HepRandomEngine * rndmEngine, const EventContext& ctx);
92
93 // some diagnostics methods are needed here too
94 float DriftTime(float zhit, const InDetDD::SiDetectorElement* element, const EventContext& ctx) const;
95 float DiffusionSigma(float zhit, const InDetDD::SiDetectorElement* element, const EventContext& ctx) const;
96 float SurfaceDriftTime(float ysurf) const;
97 float MaxDriftTime(const InDetDD::SiDetectorElement* element, const EventContext& ctx) const;
98 float MaxDiffusionSigma(const InDetDD::SiDetectorElement* element, const EventContext& ctx) const;
99
100 // methods for Taka Kondos's new charge drift model
101 // Non-const methods are used in initialization
102 void initTransportModel();
103 void init_mud_e(double T);
104 void init_mud_h(double T);
105 void initExEyArray();
106 void initPotentialValue();
107
108 void EField(double x, double y, double& Ex, double& Ey) const;
109 double induced (int istrip, double x, double y) const;
110 bool electron(double x_e, double y_e, double& vx_e, double& vy_e, double& D_e) const;
111 bool hole( double x_h, double y_h, double& vx_h, double& vy_h, double& D_h);
112 double mud_e(double E) const;
113 double mud_h(double E) const;
114 static double ExValue150(int ix, int iy) ;
115 static double EyValue150(int ix, int iy) ;
116 static double GetPotentialValue(int ix, int iy) ;
117 void holeTransport(double& x0, double& y0, double* Q_m2, double* Q_m1, double* Q_00, double* Q_p1, double* Q_p2, CLHEP::HepRandomEngine * rndmEngine);
118 void electronTransport(double& x0, double& y0, double* Q_m2, double* Q_m1, double* Q_00, double* Q_p1, double* Q_p2, CLHEP::HepRandomEngine * rndmEngine) const;
119 double inducedCharge(int& istrip, double& x, double& y, double& t) const;
120
121 IntegerProperty m_numberOfCharges{this, "NumberOfCharges", 1};
122 FloatProperty m_smallStepLength{this, "SmallStepLength", 5.};
123
125 FloatProperty m_tSurfaceDrift{this, "SurfaceDriftTime", 10., "max surface drift time"};
126
127 FloatProperty m_tfix{this, "FixedTime", -999., "ixed time"};
128 FloatProperty m_tsubtract{this, "SubtractTime", -999., "subtract drift time from mid gap"};
129
130 // -- Charge Trapping -- //
131 BooleanProperty m_doHistoTrap{this, "doHistoTrap", false, "Allow filling of histos for charge trapping effect"};
132 BooleanProperty m_doTrapping{this, "doTrapping", false, "Simulation of charge trapping effect"};
133 DoubleProperty m_Fluence{this, "Fluence", 0., "Fluence for charge trapping effect"};
134
135 //TK model settings
136 IntegerProperty m_chargeDriftModel{this, "ChargeDriftModel", ehTransport, "use enum ChargeDriftMode"};
137 IntegerProperty m_eFieldModel{this, "EFieldModel", FEMsolution, "use enum EFieldModel"};
138
139 //------TK parameters given externally by jobOptions ------------------
140 DoubleProperty m_depletionVoltage{this, "DepletionVoltage", 70.};
141 DoubleProperty m_biasVoltage{this, "BiasVoltage", 150.};
142 DoubleProperty m_magneticField{this, "MagneticField", -2.0};
143 DoubleProperty m_sensorTemperature{this, "SensorTemperature", 0.+CLHEP::STP_Temperature};
144 DoubleProperty m_transportTimeStep{this, "TransportTimeStep", 0.25};
145 DoubleProperty m_transportTimeMax{this, "TransportTimeMax", 25.0};
146
147 BooleanProperty m_isOverlay{this, "isOverlay", false, "flag for overlay"};
148
149 //ToolHandles
150 ToolHandle<ISiPropertiesTool> m_siPropertiesTool{this, "SiPropertiesTool", "SCT_SiPropertiesTool", "Tool to retrieve SCT silicon properties"};
151 ToolHandle<ISiliconConditionsTool> m_siConditionsTool{this, "SiConditionsTool", "SCT_SiliconConditionsTool", "Tool to retrieve SCT silicon information"};
152 ToolHandle<ISiLorentzAngleTool> m_lorentzAngleTool{this, "LorentzAngleTool", "SiLorentzAngleTool/SCTLorentzAngleTool", "Tool to retreive Lorentz angle"};
153
154 ServiceHandle<ITHistSvc> m_thistSvc{this, "THistSvc", "THistSvc"};
155
156 float m_tHalfwayDrift{0.};
157 float m_distInterStrip{1.0};
159
160 bool m_SurfaceDriftFlag{false};
161
162 //------TK parameters mostly fixed but can be changed externally ------------
163 double m_bulk_depth{0.0285}; //<!285 micron, expressed in cm units
164 double m_strip_pitch{0.0080}; //<! 80 micron, expressed in cm units
165 double m_depletion_depth{0.0285};
166 double m_y_origin_min{0.0}; //<! zero unless under-depleted
167
168 //------TK parameters for e, h transport --------------------------------
169 double m_kB{1.38E-23}; //<! Boltzmann const [m^2*kg/s^2/K]
170 double m_e{1.602E-19}; //<! electron charge [Coulomb]
171 double m_vs_e{11615084.7393}; //<! mobility at 273.15K
172 double m_Ec_e{6034.20429};
173 double m_vs_h{8761659.83530}; //<! hole mobility at 273.15K
174 double m_Ec_h{15366.52650};
175 double m_beta_e{};
176 double m_beta_h{};
178
179 //------TK arrays of FEM analysis -----------------------------------
180 double m_PotentialValue[81][115]{{0.}};
181 double m_ExValue150[17][115]{{0.}};
182 double m_EyValue150[17][115]{{0.}};
183
184 //------TK parameters for charge map, uses file storage of map....
185 // This member makes the class very large --- large enough that it fails
186 // ubsan's sanity checks and produces a false positive. However, it is not
187 // actually used, so comment it out. If it is ever actually needed,
188 // then it should be allocated dynamically rather than being allocated
189 // inline to the class.
190 // double m_stripCharge[5][81][285][50]{{{{0.}}}};
191 // sroe: the following were never initialised before, which begs the question:
192 // Did this code *ever* work? Has it *ever* been used?
196
197 // -- Histograms
198 TProfile* m_h_efieldz{nullptr};
199 TProfile2D* m_h_yzRamo{nullptr};
200 TProfile2D* m_h_yzEfield{nullptr};
201 TProfile2D* m_h_yEfield{nullptr};
202 TProfile2D* m_h_zEfield{nullptr};
203};
204
205#endif // SCT_DETAILEDSURFACECHARGESGENERATOR_H
#define y
#define x
Class to hold geometrical description of a silicon detector element.
ToolHandle< ISiliconConditionsTool > m_siConditionsTool
bool electron(double x_e, double y_e, double &vx_e, double &vy_e, double &D_e) const
FloatProperty m_tSurfaceDrift
related to the surface drift
virtual StatusCode finalize() override
AlgTool finalize.
virtual StatusCode initialize() override
AlgTool initialize.
virtual ~SCT_DetailedSurfaceChargesGenerator()=default
Destructor.
float MaxDriftTime(const InDetDD::SiDetectorElement *element, const EventContext &ctx) const
max drift charge equivalent to the detector thickness
float m_distInterStrip
Inter strip distance normalized to 1.
float MaxDiffusionSigma(const InDetDD::SiDetectorElement *element, const EventContext &ctx) const
max sigma diffusion
void processSiHit(const InDetDD::SiDetectorElement *element, const SiHit &phit, ISiSurfaceChargesInserter &inserter, const float eventTime, const unsigned short eventID, CLHEP::HepRandomEngine *rndmEngine, const EventContext &ctx)
SCT_DetailedSurfaceChargesGenerator(const std::string &type, const std::string &name, const IInterface *parent)
constructor
void EField(double x, double y, double &Ex, double &Ey) const
virtual void setFixedTime(float fixedTime) override
float m_distHalfInterStrip
Half way distance inter strip.
void electronTransport(double &x0, double &y0, double *Q_m2, double *Q_m1, double *Q_00, double *Q_p1, double *Q_p2, CLHEP::HepRandomEngine *rndmEngine) const
FloatProperty m_smallStepLength
max internal step along the larger G4 step
void holeTransport(double &x0, double &y0, double *Q_m2, double *Q_m1, double *Q_00, double *Q_p1, double *Q_p2, CLHEP::HepRandomEngine *rndmEngine)
double induced(int istrip, double x, double y) const
float DiffusionSigma(float zhit, const InDetDD::SiDetectorElement *element, const EventContext &ctx) const
calculate diffusion sigma from a gaussian dist scattered charge
float SurfaceDriftTime(float ysurf) const
Calculate of the surface drift time.
float DriftTime(float zhit, const InDetDD::SiDetectorElement *element, const EventContext &ctx) const
calculate drift time perpandicular to the surface for a charge at distance zhit from mid gap
static double ExValue150(int ix, int iy)
bool hole(double x_h, double y_h, double &vx_h, double &vy_h, double &D_h)
static double EyValue150(int ix, int iy)
double inducedCharge(int &istrip, double &x, double &y, double &t) const
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
const std::string process
Message Stream Member.