CMS 3D CMS Logo

SiLinearChargeDivider.cc
Go to the documentation of this file.
6 #include <fstream>
7 #include <string>
8 
10  : // Run APV in peak instead of deconvolution mode, which degrades the time resolution.
11  peakMode(conf.getParameter<bool>("APVpeakmode")),
12  // Enable interstrip Landau fluctuations within a cluster.
13  fluctuateCharge(conf.getParameter<bool>("LandauFluctuations")),
14  // Number of segments per strip into which charge is divided during
15  // simulation. If large, precision of simulation improves.
16  chargedivisionsPerStrip(conf.getParameter<int>("chargeDivisionsPerStrip")),
17  // delta cutoff in MeV, has to be same as in Geant (0.120425 MeV corresponding to 100um range for electrons)
18  deltaCut(conf.getParameter<double>("DeltaProductionCut")),
19  //Offset for digitization during the MTCC and in general for taking cosmic particle
20  //The value to be used it must be evaluated and depend on the volume defnition used
21  //for the cosimc generation (Considering only the tracker the value is 11 ns)
22  cosmicShift(conf.getUntrackedParameter<double>("CosmicDelayShift")),
23  theParticleDataTable(nullptr),
24  // Geant4 engine used to fluctuate the charge from segment to segment
25  fluctuate(new SiG4UniversalFluctuation())
26 
27 {
28  readPulseShape(conf.getParameter<edm::FileInPath>(peakMode ? "APVShapePeakFile" : "APVShapeDecoFile").fullPath());
29 }
30 
32 
33 void SiLinearChargeDivider::readPulseShape(const std::string& pulseShapeFileName) {
34  // Pulse shape file format: empty lines and comments (lines starting with '#') are ignored
35  // one line "resolution: value" is interpreted as the resolution
36  // all other lines are read as consecutive values of the shape
37  std::ifstream shapeFile(pulseShapeFileName.c_str());
38  if (!shapeFile.good()) {
39  throw cms::Exception("FileError") << "Problem opening APV Shape file: " << pulseShapeFileName;
40  }
41  pulseResolution = -1.;
43  const std::string resoPrefix{"resolution: "};
44  while (std::getline(shapeFile, line)) {
45  if ((!line.empty()) && (line.substr(1) != "#")) {
46  std::istringstream lStr{line};
47  if (line.substr(0, resoPrefix.size()) == resoPrefix) {
48  lStr.seekg(resoPrefix.size());
49  lStr >> pulseResolution;
50  } else {
51  double value;
52  while (lStr >> value) {
53  pulseValues.push_back(value);
54  }
55  }
56  }
57  }
58  if (pulseValues.empty() || (pulseResolution == -1.)) {
59  throw cms::Exception("WrongAPVPulseShape") << "Problem reading from APV pulse shape file " << pulseShapeFileName
60  << ": " << (pulseValues.empty() ? "no values" : "no resolution");
61  }
62  const auto maxIt = std::max_element(pulseValues.begin(), pulseValues.end());
63  if (std::abs((*maxIt) - 1.) > std::numeric_limits<double>::epsilon()) {
64  throw cms::Exception("WrongAPVPulseShape")
65  << "The max value of the APV pulse shape stored in the text file used in "
66  "SimGeneral/MixingModule/python/SiStripSimParameters_cfi.py is not equal to 1. Need to be fixed.";
67  }
68  pulset0Idx = std::distance(pulseValues.begin(), maxIt);
69 }
70 
72  const LocalVector& driftdir,
73  double moduleThickness,
74  const StripGeomDetUnit& det,
75  CLHEP::HepRandomEngine* engine) {
76  // signal after pulse shape correction
77  float const decSignal = TimeResponse(hit, det);
78 
79  // if out of time go home!
80  if (0 == decSignal)
81  return ionization_type();
82 
83  // Get the nass if the particle, in MeV.
84  // Protect from particles with Mass = 0, assuming then the pion mass
85  assert(theParticleDataTable != nullptr);
86  ParticleData const* particle = theParticleDataTable->particle(hit->particleType());
87  double const particleMass = particle ? particle->mass() * 1000 : 139.57;
88  double const particleCharge = particle ? particle->charge() : 1.;
89 
90  if (!particle) {
91  LogDebug("SiLinearChargeDivider") << "Cannot find particle of type " << hit->particleType()
92  << " in the PDT we assign to this particle the mass and charge of the Pion";
93  }
94 
95  int NumberOfSegmentation =
96  // if neutral: just one deposit....
97  (fabs(particleMass) < 1.e-6 || particleCharge == 0)
98  ? 1
99  :
100  // computes the number of segments from number of segments per strip times number of strips.
101  (int)(1 + chargedivisionsPerStrip *
102  fabs(driftXPos(hit->exitPoint(), driftdir, moduleThickness) -
103  driftXPos(hit->entryPoint(), driftdir, moduleThickness)) /
104  det.specificTopology().localPitch(hit->localPosition()));
105 
106  // Eloss in GeV
107  float eLoss = hit->energyLoss();
108 
109  // Prepare output
110  ionization_type _ionization_points;
111  _ionization_points.resize(NumberOfSegmentation);
112 
113  // Fluctuate charge in track subsegments
114  LocalVector direction = hit->exitPoint() - hit->entryPoint();
115  if (NumberOfSegmentation <= 1) {
116  // here I need a random... not 0.5
117  _ionization_points[0] = EnergyDepositUnit(eLoss * decSignal / eLoss, hit->entryPoint() + 0.5f * direction);
118  } else {
119  float eLossVector[NumberOfSegmentation];
120  if (fluctuateCharge) {
121  fluctuateEloss(particleMass, hit->pabs(), eLoss, direction.mag(), NumberOfSegmentation, eLossVector, engine);
122  // Save the energy of each segment
123  for (int i = 0; i != NumberOfSegmentation; i++) {
124  // take energy value from vector eLossVector,
125  _ionization_points[i] =
126  EnergyDepositUnit(eLossVector[i] * decSignal / eLoss,
127  hit->entryPoint() + float((i + 0.5) / NumberOfSegmentation) * direction);
128  }
129  } else {
130  // Save the energy of each segment
131  for (int i = 0; i != NumberOfSegmentation; i++) {
132  // take energy value from eLoss average over n.segments.
133  _ionization_points[i] =
134  EnergyDepositUnit(decSignal / float(NumberOfSegmentation),
135  hit->entryPoint() + float((i + 0.5) / NumberOfSegmentation) * direction);
136  }
137  }
138  }
139  return _ionization_points;
140 }
141 
142 void SiLinearChargeDivider::fluctuateEloss(double particleMass,
143  float particleMomentum,
144  float eloss,
145  float length,
146  int NumberOfSegs,
147  float elossVector[],
148  CLHEP::HepRandomEngine* engine) {
149  // Generate charge fluctuations.
150  float sum = 0.;
151  double deltaCutoff;
152  double mom = particleMomentum * 1000.;
153  double seglen = length / NumberOfSegs * 10.;
154  double segeloss = (1000. * eloss) / NumberOfSegs;
155  for (int i = 0; i < NumberOfSegs; i++) {
156  // The G4 routine needs momentum in MeV, mass in MeV, delta-cut in MeV,
157  // track segment length in mm, segment eloss in MeV
158  // Returns fluctuated eloss in MeV
159  // the cutoff is sometimes redefined inside, so fix it.
160  deltaCutoff = deltaCut;
161  sum += (elossVector[i] =
162  fluctuate->SampleFluctuations(mom, particleMass, deltaCutoff, seglen, segeloss, engine) / 1000.);
163  }
164 
165  if (sum > 0.) { // If fluctuations give eloss>0.
166  // Rescale to the same total eloss
167  float ratio = eloss / sum;
168  for (int ii = 0; ii < NumberOfSegs; ii++)
169  elossVector[ii] = ratio * elossVector[ii];
170  } else { // If fluctuations gives 0 eloss
171  float averageEloss = eloss / NumberOfSegs;
172  for (int ii = 0; ii < NumberOfSegs; ii++)
173  elossVector[ii] = averageEloss;
174  }
175  return;
176 }
177 
179  // x is difference between the tof and the tof for a photon (reference)
180  // converted into a bin number
181  const auto dTOF = det.surface().toGlobal(hit->localPosition()).mag() / 30. + cosmicShift - hit->tof();
182  const int x = int(dTOF / pulseResolution) + pulset0Idx;
183  if (x < 0 || x >= int(pulseValues.size()))
184  return 0;
185  return hit->energyLoss() * pulseValues[std::size_t(x)];
186 }
Vector3DBase< float, LocalTag >
SiLinearChargeDivider::SiLinearChargeDivider
SiLinearChargeDivider(const edm::ParameterSet &conf)
Definition: SiLinearChargeDivider.cc:9
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
SiLinearChargeDivider::theParticleDataTable
const ParticleDataTable * theParticleDataTable
Definition: SiLinearChargeDivider.h:49
SiLinearChargeDivider::fluctuate
std::unique_ptr< SiG4UniversalFluctuation > fluctuate
Definition: SiLinearChargeDivider.h:55
SiLinearChargeDivider::cosmicShift
const double cosmicShift
Definition: SiLinearChargeDivider.h:48
cms::cuda::assert
assert(be >=bs)
HLT_2018_cff.distance
distance
Definition: HLT_2018_cff.py:6417
DDAxes::x
SiLinearChargeDivider::pulset0Idx
unsigned int pulset0Idx
Definition: SiLinearChargeDivider.h:51
simKBmtfDigis_cfi.eLoss
eLoss
Definition: simKBmtfDigis_cfi.py:9
StripTopology::localPitch
virtual float localPitch(const LocalPoint &) const =0
SiLinearChargeDivider::~SiLinearChargeDivider
~SiLinearChargeDivider() override
Definition: SiLinearChargeDivider.cc:31
SiLinearChargeDivider::deltaCut
const double deltaCut
Definition: SiLinearChargeDivider.h:47
ParticleData
HepPDT::ParticleData ParticleData
Definition: ParticleDataTable.h:9
geometryDiff.epsilon
int epsilon
Definition: geometryDiff.py:26
edm::FileInPath
Definition: FileInPath.h:64
GeomDet::surface
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
SiLinearChargeDivider::divide
SiChargeDivider::ionization_type divide(const PSimHit *, const LocalVector &, double, const StripGeomDetUnit &det, CLHEP::HepRandomEngine *) override
Definition: SiLinearChargeDivider.cc:71
SiLinearChargeDivider::peakMode
const bool peakMode
Definition: SiLinearChargeDivider.h:44
Surface::toGlobal
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
SiLinearChargeDivider::fluctuateCharge
const bool fluctuateCharge
Definition: SiLinearChargeDivider.h:45
particleFlowDisplacedVertex_cfi.ratio
ratio
Definition: particleFlowDisplacedVertex_cfi.py:93
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
SiLinearChargeDivider.h
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
edm::ParameterSet
Definition: ParameterSet.h:36
SiLinearChargeDivider::chargedivisionsPerStrip
const int chargedivisionsPerStrip
Definition: SiLinearChargeDivider.h:46
createfilelist.int
int
Definition: createfilelist.py:10
value
Definition: value.py:1
SiLinearChargeDivider::driftXPos
float driftXPos(const Local3DPoint &pos, const LocalVector &drift, double thickness)
Definition: SiLinearChargeDivider.h:57
SiLinearChargeDivider::readPulseShape
void readPulseShape(const std::string &pulseShapeFileName)
Definition: SiLinearChargeDivider.cc:33
SiChargeDivider::ionization_type
std::vector< EnergyDepositUnit > ionization_type
Definition: SiChargeDivider.h:21
SiLinearChargeDivider::pulseResolution
double pulseResolution
Definition: SiLinearChargeDivider.h:50
PV3DBase::mag
T mag() const
Definition: PV3DBase.h:64
SiLinearChargeDivider::fluctuateEloss
void fluctuateEloss(double const particleMass, float momentum, float eloss, float length, int NumberOfSegmentation, float elossVector[], CLHEP::HepRandomEngine *)
Definition: SiLinearChargeDivider.cc:142
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
mag
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: Basic3DVectorLD.h:127
SiLinearChargeDivider::pulseValues
std::vector< double > pulseValues
Definition: SiLinearChargeDivider.h:52
LocalPoint.h
relativeConstraints.value
value
Definition: relativeConstraints.py:53
Exception
Definition: hltDiff.cc:246
SiG4UniversalFluctuation
Definition: SiG4UniversalFluctuation.h:25
StripGeomDetUnit::specificTopology
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
Definition: StripGeomDetUnit.cc:17
SiLinearChargeDivider::TimeResponse
float TimeResponse(const PSimHit *hit, const StripGeomDetUnit &det)
Definition: SiLinearChargeDivider.cc:178
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
PSimHit
Definition: PSimHit.h:15
mps_splice.line
line
Definition: mps_splice.py:76
LocalVector.h
cuy.ii
ii
Definition: cuy.py:590
EnergyDepositUnit
Definition: EnergyDepositUnit.h:10
edm::FileInPath::fullPath
std::string fullPath() const
Definition: FileInPath.cc:163
hit
Definition: SiStripHitEffFromCalibTree.cc:88
StripGeomDetUnit
Definition: StripGeomDetUnit.h:15