10 #include "CLHEP/Random/RandGaussQ.h"
42 constexpr
double m_pion = 139.571;
43 constexpr
double m_kaon = 493.677;
44 constexpr
double m_electron = 0.511;
45 constexpr
double m_muon = 105.658;
46 constexpr
double m_proton = 938.272;
50 float calcQ(
float x) {
51 constexpr
float p1 = 12.5f;
52 constexpr
float p2 = 0.2733f;
53 constexpr
float p3 = 0.147f;
56 return 0.5f * (1.f - std::copysign(
std::sqrt(1.
f - unsafe_expf<4>(-
xx * (1.
f +
p2 / (1.
f +
p3 *
xx)))), x));
62 makeDigiSimLinks_(conf_common.getUntrackedParameter<
bool>(
"makeDigiSimLinks",
true)),
63 use_ineff_from_db_(conf_specific.getParameter<
bool>(
"Inefficiency_DB")),
64 use_module_killing_(conf_specific.getParameter<
bool>(
"KillModules")),
65 use_deadmodule_DB_(conf_specific.getParameter<
bool>(
"DeadModules_DB")),
67 use_LorentzAngle_DB_(conf_specific.getParameter<
bool>(
"LorentzAngle_DB")),
70 deadModules_(use_deadmodule_DB_ ?
Parameters() : conf_specific.getParameter<
Parameters>(
"DeadModules")),
74 GeVperElectron_(3.61E-09),
75 alpha2Order_(conf_specific.getParameter<
bool>(
"Alpha2Order")),
76 addXtalk_(conf_specific.getParameter<
bool>(
"AddXTalk")),
78 interstripCoupling_(conf_specific.getParameter<double>(
"InterstripCoupling")),
80 Sigma0_(conf_specific.getParameter<double>(
"SigmaZero")),
81 SigmaCoeff_(conf_specific.getParameter<double>(
"SigmaCoeff")),
87 clusterWidth_(conf_specific.getParameter<double>(
"ClusterWidth")),
93 thePhase2ReadoutMode_(conf_specific.getParameter<
int>(
"Phase2ReadoutMode")),
100 theElectronPerADC_(conf_specific.getParameter<double>(
"ElectronPerAdc")),
103 theAdcFullScale_(conf_specific.getParameter<
int>(
"AdcFullScale")),
107 theNoiseInElectrons_(conf_specific.getParameter<double>(
"NoiseInElectrons")),
110 theReadoutNoise_(conf_specific.getParameter<double>(
"ReadoutNoiseInElec")),
115 theThresholdInE_Endcap_(conf_specific.getParameter<double>(
"ThresholdInElectrons_Endcap")),
116 theThresholdInE_Barrel_(conf_specific.getParameter<double>(
"ThresholdInElectrons_Barrel")),
119 theThresholdSmearing_Endcap_(conf_specific.getParameter<double>(
"ThresholdSmearing_Endcap")),
120 theThresholdSmearing_Barrel_(conf_specific.getParameter<double>(
"ThresholdSmearing_Barrel")),
123 theHIPThresholdInE_Endcap_(conf_specific.getParameter<double>(
"HIPThresholdInElectrons_Endcap")),
124 theHIPThresholdInE_Barrel_(conf_specific.getParameter<double>(
"HIPThresholdInElectrons_Barrel")),
127 theTofLowerCut_(conf_specific.getParameter<double>(
"TofLowerCut")),
128 theTofUpperCut_(conf_specific.getParameter<double>(
"TofUpperCut")),
131 tanLorentzAnglePerTesla_Endcap_(
132 use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>(
"TanLorentzAnglePerTesla_Endcap")),
133 tanLorentzAnglePerTesla_Barrel_(
134 use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>(
"TanLorentzAnglePerTesla_Barrel")),
137 addNoise_(conf_specific.getParameter<
bool>(
"AddNoise")),
140 addNoisyPixels_(conf_specific.getParameter<
bool>(
"AddNoisyPixels")),
143 fluctuateCharge_(conf_specific.getUntrackedParameter<
bool>(
"FluctuateCharge",
true)),
146 addPixelInefficiency_(conf_specific.getParameter<
bool>(
"AddInefficiency")),
149 addThresholdSmearing_(conf_specific.getParameter<
bool>(
"AddThresholdSmearing")),
152 pseudoRadDamage_(conf_specific.exists(
"PseudoRadDamage") ? conf_specific.getParameter<double>(
"PseudoRadDamage")
154 pseudoRadDamageRadius_(conf_specific.exists(
"PseudoRadDamageRadius")
155 ? conf_specific.getParameter<double>(
"PseudoRadDamageRadius")
161 tMax_(conf_common.getParameter<double>(
"DeltaProductionCut")),
163 badPixels_(conf_specific.getParameter<
Parameters>(
"CellsToKill")),
167 theSiPixelGainCalibrationService_(
169 subdetEfficiencies_(conf_specific) {
170 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
171 <<
"Phase2TrackerDigitizerAlgorithm constructed\n"
172 <<
"Configuration parameters:\n"
173 <<
"Threshold/Gain = "
176 <<
" ADC Scale (in bits) " <<
theAdcFullScale_ <<
" The delta cut-off is set to " <<
tMax_ <<
" pix-inefficiency "
181 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
"Phase2TrackerDigitizerAlgorithm deleted";
189 std::vector<PSimHit>::const_iterator inputEnd,
190 const size_t inputBeginGlobalIndex,
191 const uint32_t tofBin,
197 size_t simHitGlobalIndex = inputBeginGlobalIndex;
200 std::vector<PSimHit> matchedSimHits;
201 std::copy_if(inputBegin, inputEnd, std::back_inserter(matchedSimHits), [detId](
auto const&
hit) ->
bool {
202 return hit.detUnitId() == detId;
205 for (
auto const&
hit : matchedSimHits) {
206 LogDebug(
"Phase2DigitizerAlgorithm") <<
hit.particleType() <<
" " <<
hit.pabs() <<
" " <<
hit.energyLoss() <<
" "
207 <<
hit.tof() <<
" " <<
hit.trackId() <<
" " <<
hit.processType() <<
" "
208 <<
hit.detUnitId() <<
hit.entryPoint() <<
" " <<
hit.exitPoint();
209 double signalScale = 1.0;
215 const auto& collection_points =
drift(
hit, pixdet, bfield, ionization_points);
233 constexpr
float SegmentLength = 0.0010;
238 float length = direction.
mag();
240 int NumberOfSegments = static_cast<int>(length / SegmentLength);
241 if (NumberOfSegments < 1)
242 NumberOfSegments = 1;
243 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
244 <<
"enter primary_ionzation " << NumberOfSegments <<
" shift = " <<
hit.exitPoint().
x() -
hit.entryPoint().
x()
245 <<
" " <<
hit.exitPoint().
y() -
hit.entryPoint().
y() <<
" " <<
hit.exitPoint().
z() -
hit.entryPoint().
z() <<
" "
246 <<
hit.particleType() <<
" " <<
hit.pabs();
248 std::vector<float> elossVector;
249 elossVector.reserve(NumberOfSegments);
254 float averageEloss =
eLoss / NumberOfSegments;
255 elossVector.resize(NumberOfSegments, averageEloss);
258 std::vector<DigitizerUtility::EnergyDepositUnit> ionization_points;
259 ionization_points.reserve(NumberOfSegments);
261 for (
size_t i = 0;
i < elossVector.size(); ++
i) {
267 ionization_points.push_back(edu);
268 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
269 <<
i <<
" " << edu.
x() <<
" " << edu.
y() <<
" " << edu.
z() <<
" " << edu.
energy();
271 return ionization_points;
280 int pid,
float particleMomentum,
float eloss,
float length,
int NumberOfSegs)
const {
281 double particleMass = ::m_pion;
285 particleMass = ::m_electron;
287 particleMass = ::m_muon;
289 particleMass = ::m_kaon;
290 else if (pid == 2212)
291 particleMass = ::m_proton;
294 float segmentLength = length / NumberOfSegs;
298 double segmentEloss = (1000. * eloss) / NumberOfSegs;
299 std::vector<float> elossVector;
300 elossVector.reserve(NumberOfSegs);
301 for (
int i = 0;
i < NumberOfSegs; ++
i) {
305 double deltaCutoff =
tMax_;
306 float de =
fluctuate_->SampleFluctuations(particleMomentum * 1000.,
313 elossVector.push_back(de);
318 float ratio = eloss / sum;
320 std::begin(elossVector),
std::end(elossVector), std::begin(elossVector), [&
ratio](
auto const&
c) ->
float {
324 float averageEloss = eloss / NumberOfSegs;
325 elossVector.resize(NumberOfSegs, averageEloss);
340 const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points)
const {
341 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
"enter drift ";
343 std::vector<DigitizerUtility::SignalPoint> collection_points;
344 collection_points.reserve(ionization_points.size());
346 if (driftDir.
z() == 0.) {
347 LogWarning(
"Phase2TrackerDigitizerAlgorithm") <<
" pxlx: drift in z is zero ";
348 return collection_points;
351 float TanLorenzAngleX = driftDir.
x();
352 float TanLorenzAngleY = 0.;
353 float dir_z = driftDir.
z();
355 float CosLorenzAngleY = 1.;
357 TanLorenzAngleY = driftDir.
y();
364 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
365 <<
" Lorentz Tan " << TanLorenzAngleX <<
" " << TanLorenzAngleY <<
" " << CosLorenzAngleX <<
" "
366 << CosLorenzAngleY <<
" " << moduleThickness * TanLorenzAngleX <<
" " << driftDir;
368 for (
auto const&
val : ionization_points) {
370 float SegX =
val.x(), SegY =
val.y(), SegZ =
val.z();
377 float driftDistance = moduleThickness / 2. - (dir_z * SegZ);
379 if (driftDistance < 0.)
381 else if (driftDistance > moduleThickness)
382 driftDistance = moduleThickness;
385 float XDriftDueToMagField = driftDistance * TanLorenzAngleX;
386 float YDriftDueToMagField = driftDistance * TanLorenzAngleY;
389 float CloudCenterX = SegX + XDriftDueToMagField;
390 float CloudCenterY = SegY + YDriftDueToMagField;
399 float Sigma =
std::sqrt(driftLength / moduleThickness) *
Sigma0_ * moduleThickness / 0.0300;
406 float Sigma_x = Sigma / CosLorenzAngleX;
407 float Sigma_y = Sigma / CosLorenzAngleY;
410 float energyOnCollector =
val.energy();
417 energyOnCollector *=
exp(-1 * kValue * driftDistance / moduleThickness);
420 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
421 <<
"Dift DistanceZ = " << driftDistance <<
" module thickness = " << moduleThickness
422 <<
" Start Energy = " <<
val.energy() <<
" Energy after loss= " << energyOnCollector;
426 collection_points.push_back(sp);
428 return collection_points;
436 const size_t hitIndex,
437 const uint32_t tofBin,
439 const std::vector<DigitizerUtility::SignalPoint>& collection_points) {
446 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
447 <<
" enter induce_signal, " << topol->
pitch().first <<
" " << topol->
pitch().second;
450 using hit_map_type = std::map<int, float, std::less<int> >;
451 hit_map_type hit_signal;
455 for (
auto const&
v : collection_points) {
456 float CloudCenterX =
v.position().x();
457 float CloudCenterY =
v.position().y();
462 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" cloud " <<
v.position().x() <<
" " <<
v.position().y() <<
" "
463 <<
v.sigma_x() <<
" " <<
v.sigma_y() <<
" " <<
v.amplitude();
484 int IPixRightUpX = static_cast<int>(std::floor(mp.
x()));
485 int IPixRightUpY = static_cast<int>(std::floor(mp.
y()));
486 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
487 <<
" right-up " << PointRightUp <<
" " << mp.
x() <<
" " << mp.
y() <<
" " << IPixRightUpX <<
" " << IPixRightUpY;
490 int IPixLeftDownX = static_cast<int>(std::floor(mp.
x()));
491 int IPixLeftDownY = static_cast<int>(std::floor(mp.
y()));
492 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" left-down " << PointLeftDown <<
" " << mp.
x() <<
" " << mp.
y()
493 <<
" " << IPixLeftDownX <<
" " << IPixLeftDownY;
497 int numRows = topol->
nrows();
499 IPixRightUpX = numRows > IPixRightUpX ? IPixRightUpX : numRows - 1;
500 IPixRightUpY = numColumns > IPixRightUpY ? IPixRightUpY : numColumns - 1;
501 IPixLeftDownX = 0 < IPixLeftDownX ? IPixLeftDownX : 0;
502 IPixLeftDownY = 0 < IPixLeftDownY ? IPixLeftDownY : 0;
506 for (
int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) {
507 float xLB, LowerBound;
510 if (ix == 0 ||
SigmaX == 0.) {
515 LowerBound = 1 - ::calcQ((xLB - CloudCenterX) /
SigmaX);
518 float xUB, UpperBound;
519 if (ix == numRows - 1 ||
SigmaX == 0.) {
524 UpperBound = 1. - ::calcQ((xUB - CloudCenterX) /
SigmaX);
526 float TotalIntegrationRange = UpperBound - LowerBound;
527 x.emplace(ix, TotalIntegrationRange);
532 for (
int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) {
533 float yLB, LowerBound;
534 if (iy == 0 ||
SigmaY == 0.) {
539 LowerBound = 1. - ::calcQ((yLB - CloudCenterY) /
SigmaY);
542 float yUB, UpperBound;
543 if (iy == numColumns - 1 ||
SigmaY == 0.) {
548 UpperBound = 1. - ::calcQ((yUB - CloudCenterY) /
SigmaY);
551 float TotalIntegrationRange = UpperBound - LowerBound;
552 y.emplace(iy, TotalIntegrationRange);
556 for (
int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) {
557 for (
int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) {
558 float ChargeFraction =
Charge *
x[ix] *
y[iy];
560 if (ChargeFraction > 0.) {
564 hit_signal[chanFired] += ChargeFraction;
571 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
572 <<
" pixel " << ix <<
" " << iy <<
" - "
573 <<
" " << chanFired <<
" " << ChargeFraction <<
" " << mp.
x() <<
" " << mp.
y() <<
" " << lp.
x() <<
" "
581 for (
auto const& hit_s : hit_signal) {
582 int chan = hit_s.first;
596 for (
auto&
s : theSignal) {
598 if ((
s.second.ampl() +
noise) < 0.)
614 int numRows = topol->
nrows();
616 for (
auto&
s : theSignal) {
617 float signalInElectrons =
s.second.ampl();
619 std::pair<int, int> hitChan;
627 s.second.set(signalInElectrons - signalInElectrons_Xtalk);
629 if (hitChan.first != 0) {
630 auto XtalkPrev = std::make_pair(hitChan.first - 1, hitChan.second);
635 if (hitChan.first < numRows - 1) {
636 auto XtalkNext = std::make_pair(hitChan.first + 1, hitChan.second);
642 for (
auto const&
l : signalNew) {
644 auto iter = theSignal.find(
chan);
645 if (iter != theSignal.end()) {
646 theSignal[
chan] +=
l.second.ampl();
664 int numRows = topol->
nrows();
666 int numberOfPixels = numRows * numColumns;
667 std::map<int, float, std::less<int> > otherPixels;
675 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
678 << otherPixels.size();
681 for (
auto const& el : otherPixels) {
682 int iy = el.first / numRows;
683 if (iy < 0 || iy > numColumns - 1)
684 LogWarning(
"Phase2TrackerDigitizerAlgorithm") <<
" error in iy " << iy;
686 int ix = el.first - iy * numRows;
687 if (ix < 0 || ix > numRows - 1)
688 LogWarning(
"Phase2TrackerDigitizerAlgorithm") <<
" error in ix " << ix;
692 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
693 <<
" Storing noise = " << el.first <<
" " << el.second <<
" " << ix <<
" " << iy <<
" " <<
chan;
695 if (theSignal[
chan] == 0)
710 float subdetEfficiency = 1.0;
715 uint32_t layerIndex = tTopo->
pxbLayer(detID);
719 uint32_t diskIndex = 2 * tTopo->
pxfDisk(detID) - tTopo->
pxfSide(detID);
724 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" enter pixel_inefficiency " << subdetEfficiency;
728 for (
auto&
s : theSignal) {
730 if (rand > subdetEfficiency) {
761 const DetId& detId)
const {
779 float alpha2 =
std::pow(lorentzAngle, 2);
781 dir_x = -(lorentzAngle * Bfield.y() + alpha2 * Bfield.z() * Bfield.x());
782 dir_y = +(lorentzAngle * Bfield.x() - alpha2 * Bfield.z() * Bfield.y());
783 dir_z = -(1 + alpha2 *
std::pow(Bfield.z(), 2));
787 float alpha2_Endcap = 0.0;
788 float alpha2_Barrel = 0.0;
797 dir_z = -(1 + alpha2_Barrel *
std::pow(Bfield.z(), 2));
803 dir_z = -(1 + alpha2_Endcap *
std::pow(Bfield.z(), 2));
809 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" The drift direction in local coordinate is " << theDriftDirection;
810 return theDriftDirection;
819 for (
auto&
s : theSignal) {
820 std::pair<int, int> ip;
841 int Dead_detID = det_m.getParameter<
int>(
"Dead_detID");
842 Module = det_m.getParameter<
std::string>(
"Module");
843 if (detid == Dead_detID) {
853 for (
auto&
s : theSignal) {
854 std::pair<int, int> ip;
860 if (Module ==
"whole" || (Module ==
"tbmA" && ip.first >= 80 && ip.first <= 159) ||
861 (Module ==
"tbmB" && ip.first <= 79))
875 for (
size_t id = 0;
id < disabledModules.size();
id++) {
876 if (detID == disabledModules[
id].DetID) {
878 badmodule = disabledModules[
id];
888 for (
auto&
s : theSignal)
893 std::vector<GlobalPixel> badrocpositions;
894 for (
size_t j = 0; j < static_cast<size_t>(ncol);
j++) {
897 for (
auto const&
p :
path) {
902 badrocpositions.push_back(global);
909 for (
auto&
s : theSignal) {
910 std::pair<int, int> ip;
916 for (
auto const&
p : badrocpositions) {
918 if (
p.row ==
k.getParameter<
int>(
"row") && ip.first ==
k.getParameter<
int>(
"row") &&
919 std::abs(ip.second -
p.col) <
k.getParameter<
int>(
"col")) {
930 auto& theSignal =
_signal[detId];
935 if (!inserted.second) {
936 throw cms::Exception(
"LogicError") <<
"Signal was already set for DetId " << detId;
942 std::map<int, DigitizerUtility::DigiSimInfo>& digi_map,
953 float theThresholdInE = 0.;
954 float theHIPThresholdInE = 0.;
989 for (
auto const&
s : theSignal) {
991 float signalInElectrons = sig_data.
ampl();
995 if (!info_list.empty())
996 hitInfo = std::max_element(info_list.begin(), info_list.end())->
second.get();
1001 info.ot_bit = signalInElectrons > theHIPThresholdInE ?
true :
false;
1004 float charge_frac =
l.first / signalInElectrons;
1006 info.simInfoList.push_back({charge_frac,
l.second.get()});
1009 digi_map.insert({
s.first,
info});
1019 const int max_limit = 10;
1031 if (temp_signal > kink_point)
1032 temp_signal = std::floor((temp_signal - kink_point) / (
pow(2, dualslope_param - 1))) + kink_point;
1035 LogTrace(
"Phase2TrackerDigitizerAlgorithm")
1036 <<
" DetId " << detID <<
" signal_in_elec " << signal_in_elec <<
" threshold " <<
threshold
1037 <<
" signal_above_thr " << signal_in_elec -
threshold <<
" temp conversion "
1039 << temp_signal <<
" signal_in_adc " << signal_in_adc;
1041 return signal_in_adc;