File indexing completed on 2023-10-25 10:05:04
0001 #include "SiHitDigitizer.h"
0002 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0003 #include "SiLinearChargeCollectionDrifter.h"
0004 #include "SiLinearChargeDivider.h"
0005
0006 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0007 #include "SiTrivialInduceChargeOnStrips.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010
0011
0012
0013 static const double CBOLTZ_over_e_SI = 1.38E-23 / 1.6E-19;
0014 static const double noDiffusionMultiplier = 1.0e-3;
0015
0016 SiHitDigitizer::SiHitDigitizer(const edm::ParameterSet& conf)
0017 : depletionVoltage(conf.getParameter<double>("DepletionVoltage")),
0018 chargeMobility(conf.getParameter<double>("ChargeMobility")),
0019 theSiChargeDivider(new SiLinearChargeDivider(conf)),
0020 theSiChargeCollectionDrifter(new SiLinearChargeCollectionDrifter(
0021 CBOLTZ_over_e_SI * chargeMobility * conf.getParameter<double>("Temperature") *
0022 (conf.getParameter<bool>("noDiffusion") ? noDiffusionMultiplier : 1.0),
0023 conf.getParameter<double>("ChargeDistributionRMS"),
0024 depletionVoltage,
0025 conf.getParameter<double>("AppliedVoltage"))),
0026 theSiInduceChargeOnStrips(new SiTrivialInduceChargeOnStrips(conf, conf.getParameter<double>("GevPerElectron"))) {}
0027
0028 SiHitDigitizer::~SiHitDigitizer() {}
0029
0030 void SiHitDigitizer::processHit(const PSimHit* hit,
0031 const StripGeomDetUnit& det,
0032 GlobalVector bfield,
0033 float langle,
0034 std::vector<float>& locAmpl,
0035 size_t& firstChannelWithSignal,
0036 size_t& lastChannelWithSignal,
0037 const TrackerTopology* tTopo,
0038 CLHEP::HepRandomEngine* engine) {
0039
0040 double moduleThickness = det.specificSurface().bounds().thickness();
0041 double timeNormalisation = (moduleThickness * moduleThickness) / (2. * depletionVoltage * chargeMobility);
0042 LocalVector driftDir = DriftDirection(&det, bfield, langle);
0043
0044
0045 theSiInduceChargeOnStrips->induce(
0046 theSiChargeCollectionDrifter->drift(theSiChargeDivider->divide(hit, driftDir, moduleThickness, det, engine),
0047 driftDir,
0048 moduleThickness,
0049 timeNormalisation),
0050 det,
0051 locAmpl,
0052 firstChannelWithSignal,
0053 lastChannelWithSignal,
0054 tTopo);
0055 }