CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiHitDigitizer.cc
Go to the documentation of this file.
1 #include "SiHitDigitizer.h"
5 
10 
11 //#define CBOLTZ (1.38E-23)
12 //#define e_SI (1.6E-19)
13 static const double CBOLTZ_over_e_SI = 1.38E-23/1.6E-19;
14 static const double noDiffusionMultiplier = 1.0e-3;
15 
16 SiHitDigitizer::SiHitDigitizer(const edm::ParameterSet& conf, CLHEP::HepRandomEngine& eng) :
17  depletionVoltage(conf.getParameter<double>("DepletionVoltage")),
18  chargeMobility(conf.getParameter<double>("ChargeMobility")),
19  theSiChargeDivider(new SiLinearChargeDivider(conf, eng)),
20  theSiChargeCollectionDrifter(new SiLinearChargeCollectionDrifter(
21  CBOLTZ_over_e_SI * chargeMobility * conf.getParameter<double>("Temperature") * (conf.getParameter<bool>("noDiffusion") ? noDiffusionMultiplier : 1.0),
22  conf.getParameter<double>("ChargeDistributionRMS"),
23  depletionVoltage,
24  conf.getParameter<double>("AppliedVoltage"))),
25  theSiInduceChargeOnStrips(new SiTrivialInduceChargeOnStrips(conf, conf.getParameter<double>("GevPerElectron"))) {
26 }
27 
29 }
30 
31 void
32 SiHitDigitizer::processHit(const PSimHit* hit, const StripGeomDetUnit& det, GlobalVector bfield,float langle,
33  std::vector<float>& locAmpl, size_t& firstChannelWithSignal, size_t& lastChannelWithSignal,
34  const TrackerTopology *tTopo){
35 
36  // Compute the drift direction for this det
37  double moduleThickness = det.specificSurface().bounds().thickness(); // active detector thicness
38  double timeNormalisation = (moduleThickness*moduleThickness)/(2.*depletionVoltage*chargeMobility);
39  LocalVector driftDir = DriftDirection(&det,bfield,langle);
40 
41  // Fully process one SimHit
44  theSiChargeDivider->divide(hit, driftDir, moduleThickness, det),
45  driftDir,moduleThickness,timeNormalisation),
46  det,locAmpl,firstChannelWithSignal,lastChannelWithSignal,tTopo);
47 }
void processHit(const PSimHit *, const StripGeomDetUnit &, GlobalVector, float, std::vector< float > &, size_t &, size_t &, const TrackerTopology *tTopo)
LocalVector DriftDirection(const StripGeomDetUnit *_detp, GlobalVector _bfield, float langle)
const double chargeMobility
const double depletionVoltage
const Bounds & bounds() const
Definition: Surface.h:128
virtual float thickness() const =0
static const double noDiffusionMultiplier
tuple conf
Definition: dbtoconf.py:185
static const double CBOLTZ_over_e_SI
std::unique_ptr< const SiInduceChargeOnStrips > theSiInduceChargeOnStrips
std::unique_ptr< SiChargeCollectionDrifter > theSiChargeCollectionDrifter
SiHitDigitizer(const edm::ParameterSet &conf, CLHEP::HepRandomEngine &)
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:38
std::unique_ptr< SiChargeDivider > theSiChargeDivider