CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
RPixLinearChargeDivider Class Reference

#include <RPixLinearChargeDivider.h>

Public Member Functions

std::vector< RPixEnergyDepositUnitdivide (const PSimHit &hit)
 
 RPixLinearChargeDivider (const edm::ParameterSet &params, CLHEP::HepRandomEngine &eng, uint32_t det_id)
 
 ~RPixLinearChargeDivider ()
 

Private Member Functions

void FluctuateEloss (int pid, double particleMomentum, double eloss, double length, int NumberOfSegs, std::vector< RPixEnergyDepositUnit > &elossVector)
 

Private Attributes

int chargedivisions_
 
double deltaCut_
 
uint32_t det_id_
 
SiG4UniversalFluctuationfluctuate
 
bool fluctuateCharge_
 
double pitch_
 
CLHEP::HepRandomEngine & rndEngine_
 
std::vector< RPixEnergyDepositUnitthe_energy_path_distribution_
 
double thickness_
 
int verbosity_
 

Detailed Description

Definition at line 13 of file RPixLinearChargeDivider.h.

Constructor & Destructor Documentation

◆ RPixLinearChargeDivider()

RPixLinearChargeDivider::RPixLinearChargeDivider ( const edm::ParameterSet params,
CLHEP::HepRandomEngine &  eng,
uint32_t  det_id 
)

Definition at line 6 of file RPixLinearChargeDivider.cc.

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 }

References chargedivisions_, deltaCut_, fluctuate, fluctuateCharge_, CalibrationSummaryClient_cfi::params, and verbosity_.

◆ ~RPixLinearChargeDivider()

RPixLinearChargeDivider::~RPixLinearChargeDivider ( )

Definition at line 17 of file RPixLinearChargeDivider.cc.

17 { delete fluctuate; }

References fluctuate.

Member Function Documentation

◆ divide()

std::vector< RPixEnergyDepositUnit > RPixLinearChargeDivider::divide ( const PSimHit hit)

Definition at line 19 of file RPixLinearChargeDivider.cc.

19  {
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("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("RPixLinearChargeDivider") << "energy dep. sum=" << sum;
57  }
58 
60 }

References chargedivisions_, det_id_, simKBmtfDigis_cfi::eLoss, fluctuateCharge_, FluctuateEloss(), mps_fire::i, PV3DBase< T, PVType, FrameType >::mag(), the_energy_path_distribution_, verbosity_, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

◆ FluctuateEloss()

void RPixLinearChargeDivider::FluctuateEloss ( int  pid,
double  particleMomentum,
double  eloss,
double  length,
int  NumberOfSegs,
std::vector< RPixEnergyDepositUnit > &  elossVector 
)
private

Definition at line 62 of file RPixLinearChargeDivider.cc.

67  {
68  double particleMass = 139.6; // Mass in MeV, Assume pion
69  pid = std::abs(pid);
70  if (pid != 211) { // Mass in MeV
71  if (pid == 11)
72  particleMass = 0.511;
73  else if (pid == 13)
74  particleMass = 105.7;
75  else if (pid == 321)
76  particleMass = 493.7;
77  else if (pid == 2212)
78  particleMass = 938.3;
79  }
80 
81  double segmentLength = length / NumberOfSegs;
82 
83  // Generate charge fluctuations.
84  double de = 0.;
85  double sum = 0.;
86  double segmentEloss = (eloss * 1000) / NumberOfSegs; //eloss in MeV
87  for (int i = 0; i < NumberOfSegs; i++) {
88  double deltaCutoff = deltaCut_;
90  particleMomentum * 1000, particleMass, deltaCutoff, segmentLength, segmentEloss, &(rndEngine_)) /
91  1000; //convert to GeV
92  elossVector[i].setEnergy(de);
93  sum += de;
94  }
95 
96  if (sum > 0.) { // If fluctuations give eloss>0.
97  // Rescale to the same total eloss
98  double ratio = eloss / sum;
99  for (int ii = 0; ii < NumberOfSegs; ii++)
100  elossVector[ii].setEnergy(ratio * elossVector[ii].Energy());
101  } else { // If fluctuations gives 0 eloss
102  double averageEloss = eloss / NumberOfSegs;
103  for (int ii = 0; ii < NumberOfSegs; ii++)
104  elossVector[ii].setEnergy(averageEloss);
105  }
106  return;
107 }

References funct::abs(), deltaCut_, EcalCondDBWriter_cfi::Energy, fluctuate, mps_fire::i, cuy::ii, particleFlowDisplacedVertex_cfi::ratio, rndEngine_, and SiG4UniversalFluctuation::SampleFluctuations().

Referenced by divide().

Member Data Documentation

◆ chargedivisions_

int RPixLinearChargeDivider::chargedivisions_
private

Definition at line 24 of file RPixLinearChargeDivider.h.

Referenced by divide(), and RPixLinearChargeDivider().

◆ deltaCut_

double RPixLinearChargeDivider::deltaCut_
private

Definition at line 25 of file RPixLinearChargeDivider.h.

Referenced by FluctuateEloss(), and RPixLinearChargeDivider().

◆ det_id_

uint32_t RPixLinearChargeDivider::det_id_
private

Definition at line 22 of file RPixLinearChargeDivider.h.

Referenced by divide().

◆ fluctuate

SiG4UniversalFluctuation* RPixLinearChargeDivider::fluctuate
private

◆ fluctuateCharge_

bool RPixLinearChargeDivider::fluctuateCharge_
private

Definition at line 23 of file RPixLinearChargeDivider.h.

Referenced by divide(), and RPixLinearChargeDivider().

◆ pitch_

double RPixLinearChargeDivider::pitch_
private

Definition at line 26 of file RPixLinearChargeDivider.h.

◆ rndEngine_

CLHEP::HepRandomEngine& RPixLinearChargeDivider::rndEngine_
private

Definition at line 21 of file RPixLinearChargeDivider.h.

Referenced by FluctuateEloss().

◆ the_energy_path_distribution_

std::vector<RPixEnergyDepositUnit> RPixLinearChargeDivider::the_energy_path_distribution_
private

Definition at line 29 of file RPixLinearChargeDivider.h.

Referenced by divide().

◆ thickness_

double RPixLinearChargeDivider::thickness_
private

Definition at line 27 of file RPixLinearChargeDivider.h.

◆ verbosity_

int RPixLinearChargeDivider::verbosity_
private

Definition at line 31 of file RPixLinearChargeDivider.h.

Referenced by divide(), and RPixLinearChargeDivider().

Vector3DBase< float, LocalTag >
mps_fire.i
i
Definition: mps_fire.py:428
RPixLinearChargeDivider::deltaCut_
double deltaCut_
Definition: RPixLinearChargeDivider.h:25
EcalCondDBWriter_cfi.Energy
Energy
Definition: EcalCondDBWriter_cfi.py:152
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
simKBmtfDigis_cfi.eLoss
eLoss
Definition: simKBmtfDigis_cfi.py:9
RPixLinearChargeDivider::fluctuate
SiG4UniversalFluctuation * fluctuate
Definition: RPixLinearChargeDivider.h:30
RPixLinearChargeDivider::FluctuateEloss
void FluctuateEloss(int pid, double particleMomentum, double eloss, double length, int NumberOfSegs, std::vector< RPixEnergyDepositUnit > &elossVector)
Definition: RPixLinearChargeDivider.cc:62
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
RPixLinearChargeDivider::fluctuateCharge_
bool fluctuateCharge_
Definition: RPixLinearChargeDivider.h:23
RPixLinearChargeDivider::chargedivisions_
int chargedivisions_
Definition: RPixLinearChargeDivider.h:24
RPixLinearChargeDivider::det_id_
uint32_t det_id_
Definition: RPixLinearChargeDivider.h:22
SiG4UniversalFluctuation::SampleFluctuations
double SampleFluctuations(const double momentum, const double mass, double &tmax, const double length, const double meanLoss, CLHEP::HepRandomEngine *)
Definition: SiG4UniversalFluctuation.cc:53
particleFlowDisplacedVertex_cfi.ratio
ratio
Definition: particleFlowDisplacedVertex_cfi.py:93
RPixLinearChargeDivider::rndEngine_
CLHEP::HepRandomEngine & rndEngine_
Definition: RPixLinearChargeDivider.h:21
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
PV3DBase::mag
T mag() const
Definition: PV3DBase.h:64
SiG4UniversalFluctuation
Definition: SiG4UniversalFluctuation.h:25
RPixLinearChargeDivider::the_energy_path_distribution_
std::vector< RPixEnergyDepositUnit > the_energy_path_distribution_
Definition: RPixLinearChargeDivider.h:29
RPixLinearChargeDivider::verbosity_
int verbosity_
Definition: RPixLinearChargeDivider.h:31
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
cuy.ii
ii
Definition: cuy.py:590
hit
Definition: SiStripHitEffFromCalibTree.cc:88