CMS 3D CMS Logo

RPixLinearChargeDivider.cc
Go to the documentation of this file.
5 
7  CLHEP::HepRandomEngine& eng,
8  uint32_t det_id)
9  : rndEngine_(eng), det_id_(det_id) {
10  verbosity_ = params.getParameter<int>("RPixVerbosity");
12  fluctuateCharge_ = params.getParameter<bool>("RPixLandauFluctuations");
13  chargedivisions_ = params.getParameter<int>("RPixChargeDivisions");
14  deltaCut_ = params.getParameter<double>("RPixDeltaProductionCut");
15 }
16 
18 
19 std::vector<RPixEnergyDepositUnit> RPixLinearChargeDivider::divide(const PSimHit& hit) {
20  LocalVector direction = hit.exitPoint() - hit.entryPoint();
21  if (direction.z() > 10 || direction.x() > 200 || direction.y() > 200) {
24  }
25 
26  int NumberOfSegmentation = chargedivisions_;
27  double eLoss = hit.energyLoss(); // Eloss in GeV
28  the_energy_path_distribution_.resize(NumberOfSegmentation);
29 
30  if (fluctuateCharge_) {
31  int pid = hit.particleType();
32  double momentum = hit.pabs();
33  double length = direction.mag(); // Track length in Silicon
34  FluctuateEloss(pid, momentum, eLoss, length, NumberOfSegmentation, the_energy_path_distribution_);
35  for (int i = 0; i < NumberOfSegmentation; i++) {
36  the_energy_path_distribution_[i].setPosition(hit.entryPoint() +
37  double((i + 0.5) / NumberOfSegmentation) * direction);
38  }
39  } else {
40  for (int i = 0; i < NumberOfSegmentation; i++) {
41  the_energy_path_distribution_[i].setPosition(hit.entryPoint() +
42  double((i + 0.5) / NumberOfSegmentation) * direction);
43  the_energy_path_distribution_[i].setEnergy(eLoss / (double)NumberOfSegmentation);
44  }
45  }
46 
47  if (verbosity_) {
48  edm::LogInfo("PPS") << "RPixLinearChargeDivider " << det_id_ << " charge along the track:";
49  double sum = 0;
50  for (unsigned int i = 0; i < the_energy_path_distribution_.size(); i++) {
51  edm::LogInfo("RPixLinearChargeDivider")
52  << the_energy_path_distribution_[i].Position().x() << " " << the_energy_path_distribution_[i].Position().y()
53  << " " << the_energy_path_distribution_[i].Position().z() << " " << the_energy_path_distribution_[i].Energy();
54  sum += the_energy_path_distribution_[i].Energy();
55  }
56  edm::LogInfo("PPS") << "RPixLinearChargeDivider "
57  << "energy dep. sum=" << sum;
58  }
59 
61 }
62 
64  double particleMomentum,
65  double eloss,
66  double length,
67  int NumberOfSegs,
68  std::vector<RPixEnergyDepositUnit>& elossVector) {
69  double particleMass = 139.6; // Mass in MeV, Assume pion
70  pid = std::abs(pid);
71  if (pid != 211) { // Mass in MeV
72  if (pid == 11)
73  particleMass = 0.511;
74  else if (pid == 13)
75  particleMass = 105.7;
76  else if (pid == 321)
77  particleMass = 493.7;
78  else if (pid == 2212)
79  particleMass = 938.3;
80  }
81 
82  double segmentLength = length / NumberOfSegs;
83 
84  // Generate charge fluctuations.
85  double de = 0.;
86  double sum = 0.;
87  double segmentEloss = (eloss * 1000) / NumberOfSegs; //eloss in MeV
88  for (int i = 0; i < NumberOfSegs; i++) {
89  double deltaCutoff = deltaCut_;
91  particleMomentum * 1000, particleMass, deltaCutoff, segmentLength, segmentEloss, &(rndEngine_)) /
92  1000; //convert to GeV
93  elossVector[i].setEnergy(de);
94  sum += de;
95  }
96 
97  if (sum > 0.) { // If fluctuations give eloss>0.
98  // Rescale to the same total eloss
99  double ratio = eloss / sum;
100  for (int ii = 0; ii < NumberOfSegs; ii++)
101  elossVector[ii].setEnergy(ratio * elossVector[ii].Energy());
102  } else { // If fluctuations gives 0 eloss
103  double averageEloss = eloss / NumberOfSegs;
104  for (int ii = 0; ii < NumberOfSegs; ii++)
105  elossVector[ii].setEnergy(averageEloss);
106  }
107  return;
108 }
RPixLinearChargeDivider(const edm::ParameterSet &params, CLHEP::HepRandomEngine &eng, uint32_t det_id)
std::vector< RPixEnergyDepositUnit > the_energy_path_distribution_
T z() const
Definition: PV3DBase.h:61
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
SiG4UniversalFluctuation * fluctuate
void FluctuateEloss(int pid, double particleMomentum, double eloss, double length, int NumberOfSegs, std::vector< RPixEnergyDepositUnit > &elossVector)
CLHEP::HepRandomEngine & rndEngine_
T mag() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ii
Definition: cuy.py:589
Log< level::Info, false > LogInfo
double SampleFluctuations(const double momentum, const double mass, double &tmax, const double length, const double meanLoss, CLHEP::HepRandomEngine *)
std::vector< RPixEnergyDepositUnit > divide(const PSimHit &hit)