CMS 3D CMS Logo

SiLinearChargeDivider.h
Go to the documentation of this file.
1 #ifndef Tracker_SiLinearChargeDivider_H
2 #define Tracker_SiLinearChargeDivider_H
3 
4 #include <memory>
5 
7 
8 #include "SiChargeDivider.h"
13 
15 
16 namespace CLHEP {
17  class HepRandomEngine;
18 }
19 
28 public:
29  // constructor
31 
32  // destructor
33  ~SiLinearChargeDivider() override;
34 
35  // main method: divide the charge (from the PSimHit) into several energy deposits in the bulk
37  const PSimHit*, const LocalVector&, double, const StripGeomDetUnit& det, CLHEP::HepRandomEngine*) override;
38 
39  // set the ParticleDataTable (used to fluctuate the charge properly)
40  void setParticleDataTable(const ParticleDataTable* pdt) override { theParticleDataTable = pdt; }
41 
42 private:
43  // configuration data
44  const bool peakMode;
45  const bool fluctuateCharge;
47  const double deltaCut;
48  const double cosmicShift;
51  unsigned int pulset0Idx;
52  std::vector<double> pulseValues;
53 
54  // Geant4 engine used by fluctuateEloss()
55  std::unique_ptr<SiG4UniversalFluctuation> fluctuate;
56  // utility: drifts the charge to the surface to estimate the number of relevant strips
57  inline float driftXPos(const Local3DPoint& pos, const LocalVector& drift, double thickness) {
58  return pos.x() + (thickness / 2. - pos.z()) * drift.x() / drift.z();
59  }
60  // fluctuate the Eloss
61  void fluctuateEloss(double const particleMass,
62  float momentum,
63  float eloss,
64  float length,
65  int NumberOfSegmentation,
66  float elossVector[],
67  CLHEP::HepRandomEngine*);
68  // time response (from the pulse shape)
69  float TimeResponse(const PSimHit* hit, const StripGeomDetUnit& det);
70  void readPulseShape(const std::string& pulseShapeFileName);
71 };
72 
73 #endif
const ParticleDataTable * theParticleDataTable
HepPDT::ParticleDataTable ParticleDataTable
T z() const
Definition: PV3DBase.h:61
void setParticleDataTable(const ParticleDataTable *pdt) override
LocalVector drift(const StripGeomDetUnit *, const MagneticField &, const SiStripLorentzAngle &)
Definition: ShallowTools.cc:36
float TimeResponse(const PSimHit *hit, const StripGeomDetUnit &det)
SiLinearChargeDivider(const edm::ParameterSet &conf)
float driftXPos(const Local3DPoint &pos, const LocalVector &drift, double thickness)
void fluctuateEloss(double const particleMass, float momentum, float eloss, float length, int NumberOfSegmentation, float elossVector[], CLHEP::HepRandomEngine *)
T x() const
Definition: PV3DBase.h:59
std::unique_ptr< SiG4UniversalFluctuation > fluctuate
void readPulseShape(const std::string &pulseShapeFileName)
std::vector< EnergyDepositUnit > ionization_type
SiChargeDivider::ionization_type divide(const PSimHit *, const LocalVector &, double, const StripGeomDetUnit &det, CLHEP::HepRandomEngine *) override
std::vector< double > pulseValues