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));
63 makeDigiSimLinks_(conf_common.getUntrackedParameter<
bool>(
"makeDigiSimLinks",
true)),
64 use_ineff_from_db_(conf_specific.getParameter<
bool>(
"Inefficiency_DB")),
65 use_module_killing_(conf_specific.getParameter<
bool>(
"KillModules")),
66 use_deadmodule_DB_(conf_specific.getParameter<
bool>(
"DeadModules_DB")),
68 use_LorentzAngle_DB_(conf_specific.getParameter<
bool>(
"LorentzAngle_DB")),
71 deadModules_(use_deadmodule_DB_ ?
Parameters() : conf_specific.getParameter<
Parameters>(
"DeadModules")),
75 GeVperElectron_(3.61E-09),
76 alpha2Order_(conf_specific.getParameter<
bool>(
"Alpha2Order")),
77 addXtalk_(conf_specific.getParameter<
bool>(
"AddXTalk")),
79 interstripCoupling_(conf_specific.getParameter<double>(
"InterstripCoupling")),
81 Sigma0_(conf_specific.getParameter<double>(
"SigmaZero")),
82 SigmaCoeff_(conf_specific.getParameter<double>(
"SigmaCoeff")),
88 clusterWidth_(conf_specific.getParameter<double>(
"ClusterWidth")),
94 thePhase2ReadoutMode_(conf_specific.getParameter<
int>(
"Phase2ReadoutMode")),
101 theElectronPerADC_(conf_specific.getParameter<double>(
"ElectronPerAdc")),
104 theAdcFullScale_(conf_specific.getParameter<
int>(
"AdcFullScale")),
108 theNoiseInElectrons_(conf_specific.getParameter<double>(
"NoiseInElectrons")),
111 theReadoutNoise_(conf_specific.getParameter<double>(
"ReadoutNoiseInElec")),
116 theThresholdInE_Endcap_(conf_specific.getParameter<double>(
"ThresholdInElectrons_Endcap")),
117 theThresholdInE_Barrel_(conf_specific.getParameter<double>(
"ThresholdInElectrons_Barrel")),
120 theThresholdSmearing_Endcap_(conf_specific.getParameter<double>(
"ThresholdSmearing_Endcap")),
121 theThresholdSmearing_Barrel_(conf_specific.getParameter<double>(
"ThresholdSmearing_Barrel")),
124 theHIPThresholdInE_Endcap_(conf_specific.getParameter<double>(
"HIPThresholdInElectrons_Endcap")),
125 theHIPThresholdInE_Barrel_(conf_specific.getParameter<double>(
"HIPThresholdInElectrons_Barrel")),
128 theTofLowerCut_(conf_specific.getParameter<double>(
"TofLowerCut")),
129 theTofUpperCut_(conf_specific.getParameter<double>(
"TofUpperCut")),
132 tanLorentzAnglePerTesla_Endcap_(
133 use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>(
"TanLorentzAnglePerTesla_Endcap")),
134 tanLorentzAnglePerTesla_Barrel_(
135 use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>(
"TanLorentzAnglePerTesla_Barrel")),
138 addNoise_(conf_specific.getParameter<
bool>(
"AddNoise")),
141 addNoisyPixels_(conf_specific.getParameter<
bool>(
"AddNoisyPixels")),
144 fluctuateCharge_(conf_specific.getUntrackedParameter<
bool>(
"FluctuateCharge",
true)),
147 addPixelInefficiency_(conf_specific.getParameter<
bool>(
"AddInefficiency")),
150 addThresholdSmearing_(conf_specific.getParameter<
bool>(
"AddThresholdSmearing")),
153 pseudoRadDamage_(conf_specific.exists(
"PseudoRadDamage") ? conf_specific.getParameter<double>(
"PseudoRadDamage")
155 pseudoRadDamageRadius_(conf_specific.exists(
"PseudoRadDamageRadius")
156 ? conf_specific.getParameter<double>(
"PseudoRadDamageRadius")
162 tMax_(conf_common.getParameter<double>(
"DeltaProductionCut")),
164 badPixels_(conf_specific.getParameter<
Parameters>(
"CellsToKill")),
168 theSiPixelGainCalibrationService_(
170 subdetEfficiencies_(conf_specific) {
171 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
172 <<
"Phase2TrackerDigitizerAlgorithm constructed\n"
173 <<
"Configuration parameters:\n"
174 <<
"Threshold/Gain = "
177 <<
" ADC Scale (in bits) " <<
theAdcFullScale_ <<
" The delta cut-off is set to " <<
tMax_ <<
" pix-inefficiency "
182 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
"Phase2TrackerDigitizerAlgorithm deleted";
190 std::vector<PSimHit>::const_iterator inputEnd,
191 const size_t inputBeginGlobalIndex,
192 const uint32_t tofBin,
198 size_t simHitGlobalIndex = inputBeginGlobalIndex;
201 std::vector<PSimHit> matchedSimHits;
202 std::copy_if(inputBegin, inputEnd, std::back_inserter(matchedSimHits), [detId](
auto const&
hit) ->
bool {
203 return hit.detUnitId() == detId;
206 for (
auto const&
hit : matchedSimHits) {
207 LogDebug(
"Phase2DigitizerAlgorithm") <<
hit.particleType() <<
" " <<
hit.pabs() <<
" " <<
hit.energyLoss() <<
" "
208 <<
hit.tof() <<
" " <<
hit.trackId() <<
" " <<
hit.processType() <<
" "
209 <<
hit.detUnitId() <<
hit.entryPoint() <<
" " <<
hit.exitPoint();
210 double signalScale = 1.0;
216 const auto& collection_points =
drift(
hit, pixdet, bfield, ionization_points);
234 constexpr
float SegmentLength = 0.0010;
239 float length = direction.
mag();
241 int NumberOfSegments = static_cast<int>(length / SegmentLength);
242 if (NumberOfSegments < 1)
243 NumberOfSegments = 1;
244 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
245 <<
"enter primary_ionzation " << NumberOfSegments <<
" shift = " <<
hit.exitPoint().
x() -
hit.entryPoint().
x()
246 <<
" " <<
hit.exitPoint().
y() -
hit.entryPoint().
y() <<
" " <<
hit.exitPoint().
z() -
hit.entryPoint().
z() <<
" "
247 <<
hit.particleType() <<
" " <<
hit.pabs();
249 std::vector<float> elossVector;
250 elossVector.reserve(NumberOfSegments);
255 float averageEloss =
eLoss / NumberOfSegments;
256 elossVector.resize(NumberOfSegments, averageEloss);
259 std::vector<DigitizerUtility::EnergyDepositUnit> ionization_points;
260 ionization_points.reserve(NumberOfSegments);
262 for (
size_t i = 0;
i < elossVector.size(); ++
i) {
268 ionization_points.push_back(edu);
269 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
270 <<
i <<
" " << edu.
x() <<
" " << edu.
y() <<
" " << edu.
z() <<
" " << edu.
energy();
272 return ionization_points;
281 int pid,
float particleMomentum,
float eloss,
float length,
int NumberOfSegs)
const {
282 double particleMass = ::m_pion;
286 particleMass = ::m_electron;
288 particleMass = ::m_muon;
290 particleMass = ::m_kaon;
291 else if (pid == 2212)
292 particleMass = ::m_proton;
295 float segmentLength = length / NumberOfSegs;
299 double segmentEloss = (1000. * eloss) / NumberOfSegs;
300 std::vector<float> elossVector;
301 elossVector.reserve(NumberOfSegs);
302 for (
int i = 0;
i < NumberOfSegs; ++
i) {
306 double deltaCutoff =
tMax_;
307 float de =
fluctuate_->SampleFluctuations(particleMomentum * 1000.,
314 elossVector.push_back(de);
319 float ratio = eloss / sum;
321 std::begin(elossVector),
std::end(elossVector), std::begin(elossVector), [&
ratio](
auto const&
c) ->
float {
325 float averageEloss = eloss / NumberOfSegs;
326 elossVector.resize(NumberOfSegs, averageEloss);
341 const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points)
const {
342 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
"enter drift ";
344 std::vector<DigitizerUtility::SignalPoint> collection_points;
345 collection_points.reserve(ionization_points.size());
347 if (driftDir.
z() == 0.) {
348 LogWarning(
"Phase2TrackerDigitizerAlgorithm") <<
" pxlx: drift in z is zero ";
349 return collection_points;
352 float TanLorenzAngleX = driftDir.
x();
353 float TanLorenzAngleY = 0.;
354 float dir_z = driftDir.
z();
356 float CosLorenzAngleY = 1.;
358 TanLorenzAngleY = driftDir.
y();
365 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
366 <<
" Lorentz Tan " << TanLorenzAngleX <<
" " << TanLorenzAngleY <<
" " << CosLorenzAngleX <<
" "
367 << CosLorenzAngleY <<
" " << moduleThickness * TanLorenzAngleX <<
" " << driftDir;
369 for (
auto const&
val : ionization_points) {
371 float SegX =
val.x(), SegY =
val.y(), SegZ =
val.z();
378 float driftDistance = moduleThickness / 2. - (dir_z * SegZ);
380 if (driftDistance < 0.)
382 else if (driftDistance > moduleThickness)
383 driftDistance = moduleThickness;
386 float XDriftDueToMagField = driftDistance * TanLorenzAngleX;
387 float YDriftDueToMagField = driftDistance * TanLorenzAngleY;
390 float CloudCenterX = SegX + XDriftDueToMagField;
391 float CloudCenterY = SegY + YDriftDueToMagField;
400 float Sigma =
std::sqrt(driftLength / moduleThickness) *
Sigma0_ * moduleThickness / 0.0300;
407 float Sigma_x = Sigma / CosLorenzAngleX;
408 float Sigma_y = Sigma / CosLorenzAngleY;
411 float energyOnCollector =
val.energy();
418 energyOnCollector *=
exp(-1 * kValue * driftDistance / moduleThickness);
421 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
422 <<
"Dift DistanceZ = " << driftDistance <<
" module thickness = " << moduleThickness
423 <<
" Start Energy = " <<
val.energy() <<
" Energy after loss= " << energyOnCollector;
427 collection_points.push_back(sp);
429 return collection_points;
437 const size_t hitIndex,
438 const uint32_t tofBin,
440 const std::vector<DigitizerUtility::SignalPoint>& collection_points) {
447 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
448 <<
" enter induce_signal, " << topol->
pitch().first <<
" " << topol->
pitch().second;
451 using hit_map_type = std::map<int, float, std::less<int> >;
452 hit_map_type hit_signal;
456 for (
auto const&
v : collection_points) {
457 float CloudCenterX =
v.position().x();
458 float CloudCenterY =
v.position().y();
463 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" cloud " <<
v.position().x() <<
" " <<
v.position().y() <<
" "
464 <<
v.sigma_x() <<
" " <<
v.sigma_y() <<
" " <<
v.amplitude();
485 int IPixRightUpX = static_cast<int>(std::floor(mp.
x()));
486 int IPixRightUpY = static_cast<int>(std::floor(mp.
y()));
487 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
488 <<
" right-up " << PointRightUp <<
" " << mp.
x() <<
" " << mp.
y() <<
" " << IPixRightUpX <<
" " << IPixRightUpY;
491 int IPixLeftDownX = static_cast<int>(std::floor(mp.
x()));
492 int IPixLeftDownY = static_cast<int>(std::floor(mp.
y()));
493 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" left-down " << PointLeftDown <<
" " << mp.
x() <<
" " << mp.
y()
494 <<
" " << IPixLeftDownX <<
" " << IPixLeftDownY;
498 int numRows = topol->
nrows();
500 IPixRightUpX = numRows > IPixRightUpX ? IPixRightUpX : numRows - 1;
501 IPixRightUpY = numColumns > IPixRightUpY ? IPixRightUpY : numColumns - 1;
502 IPixLeftDownX = 0 < IPixLeftDownX ? IPixLeftDownX : 0;
503 IPixLeftDownY = 0 < IPixLeftDownY ? IPixLeftDownY : 0;
507 for (
int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) {
508 float xLB, LowerBound;
511 if (ix == 0 ||
SigmaX == 0.) {
516 LowerBound = 1 - ::calcQ((xLB - CloudCenterX) /
SigmaX);
519 float xUB, UpperBound;
520 if (ix == numRows - 1 ||
SigmaX == 0.) {
525 UpperBound = 1. - ::calcQ((xUB - CloudCenterX) /
SigmaX);
527 float TotalIntegrationRange = UpperBound - LowerBound;
528 x.emplace(ix, TotalIntegrationRange);
533 for (
int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) {
534 float yLB, LowerBound;
535 if (iy == 0 ||
SigmaY == 0.) {
540 LowerBound = 1. - ::calcQ((yLB - CloudCenterY) /
SigmaY);
543 float yUB, UpperBound;
544 if (iy == numColumns - 1 ||
SigmaY == 0.) {
549 UpperBound = 1. - ::calcQ((yUB - CloudCenterY) /
SigmaY);
552 float TotalIntegrationRange = UpperBound - LowerBound;
553 y.emplace(iy, TotalIntegrationRange);
557 for (
int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) {
558 for (
int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) {
559 float ChargeFraction =
Charge *
x[ix] *
y[iy];
561 if (ChargeFraction > 0.) {
565 hit_signal[chanFired] += ChargeFraction;
572 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
573 <<
" pixel " << ix <<
" " << iy <<
" - "
574 <<
" " << chanFired <<
" " << ChargeFraction <<
" " << mp.
x() <<
" " << mp.
y() <<
" " << lp.
x() <<
" "
582 for (
auto const& hit_s : hit_signal) {
583 int chan = hit_s.first;
597 for (
auto&
s : theSignal) {
599 if ((
s.second.ampl() +
noise) < 0.)
615 int numRows = topol->
nrows();
617 for (
auto&
s : theSignal) {
618 float signalInElectrons =
s.second.ampl();
620 std::pair<int, int> hitChan;
628 s.second.set(signalInElectrons - signalInElectrons_Xtalk);
630 if (hitChan.first != 0) {
631 auto XtalkPrev = std::make_pair(hitChan.first - 1, hitChan.second);
636 if (hitChan.first < numRows - 1) {
637 auto XtalkNext = std::make_pair(hitChan.first + 1, hitChan.second);
643 for (
auto const&
l : signalNew) {
645 auto iter = theSignal.find(
chan);
646 if (iter != theSignal.end()) {
647 theSignal[
chan] +=
l.second.ampl();
665 int numRows = topol->
nrows();
667 int numberOfPixels = numRows * numColumns;
668 std::map<int, float, std::less<int> > otherPixels;
676 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
679 << otherPixels.size();
682 for (
auto const& el : otherPixels) {
683 int iy = el.first / numRows;
684 if (iy < 0 || iy > numColumns - 1)
685 LogWarning(
"Phase2TrackerDigitizerAlgorithm") <<
" error in iy " << iy;
687 int ix = el.first - iy * numRows;
688 if (ix < 0 || ix > numRows - 1)
689 LogWarning(
"Phase2TrackerDigitizerAlgorithm") <<
" error in ix " << ix;
693 LogDebug(
"Phase2TrackerDigitizerAlgorithm")
694 <<
" Storing noise = " << el.first <<
" " << el.second <<
" " << ix <<
" " << iy <<
" " <<
chan;
696 if (theSignal[
chan] == 0)
711 float subdetEfficiency = 1.0;
716 uint32_t layerIndex = tTopo->
pxbLayer(detID);
720 uint32_t diskIndex = 2 * tTopo->
pxfDisk(detID) - tTopo->
pxfSide(detID);
725 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" enter pixel_inefficiency " << subdetEfficiency;
729 for (
auto&
s : theSignal) {
731 if (rand > subdetEfficiency) {
762 const DetId& detId)
const {
780 float alpha2 =
std::pow(lorentzAngle, 2);
782 dir_x = -(lorentzAngle * Bfield.y() + alpha2 * Bfield.z() * Bfield.x());
783 dir_y = +(lorentzAngle * Bfield.x() - alpha2 * Bfield.z() * Bfield.y());
784 dir_z = -(1 + alpha2 *
std::pow(Bfield.z(), 2));
788 float alpha2_Endcap = 0.0;
789 float alpha2_Barrel = 0.0;
798 dir_z = -(1 + alpha2_Barrel *
std::pow(Bfield.z(), 2));
804 dir_z = -(1 + alpha2_Endcap *
std::pow(Bfield.z(), 2));
810 LogDebug(
"Phase2TrackerDigitizerAlgorithm") <<
" The drift direction in local coordinate is " << theDriftDirection;
811 return theDriftDirection;
820 for (
auto&
s : theSignal) {
821 std::pair<int, int> ip;
842 int Dead_detID = det_m.getParameter<
int>(
"Dead_detID");
843 Module = det_m.getParameter<
std::string>(
"Module");
844 if (detid == Dead_detID) {
854 for (
auto&
s : theSignal) {
855 std::pair<int, int> ip;
861 if (Module ==
"whole" || (Module ==
"tbmA" && ip.first >= 80 && ip.first <= 159) ||
862 (Module ==
"tbmB" && ip.first <= 79))
876 for (
size_t id = 0;
id < disabledModules.size();
id++) {
877 if (detID == disabledModules[
id].DetID) {
879 badmodule = disabledModules[
id];
889 for (
auto&
s : theSignal)
894 std::vector<GlobalPixel> badrocpositions;
895 for (
size_t j = 0; j < static_cast<size_t>(ncol);
j++) {
898 for (
auto const&
p :
path) {
903 badrocpositions.push_back(global);
910 for (
auto&
s : theSignal) {
911 std::pair<int, int> ip;
917 for (
auto const&
p : badrocpositions) {
919 if (
p.row ==
k.getParameter<
int>(
"row") && ip.first ==
k.getParameter<
int>(
"row") &&
920 std::abs(ip.second -
p.col) <
k.getParameter<
int>(
"col")) {
931 auto& theSignal =
_signal[detId];
936 if (!inserted.second) {
937 throw cms::Exception(
"LogicError") <<
"Signal was already set for DetId " << detId;
943 std::map<int, DigitizerUtility::DigiSimInfo>& digi_map,
954 float theThresholdInE = 0.;
955 float theHIPThresholdInE = 0.;
990 for (
auto const&
s : theSignal) {
992 float signalInElectrons = sig_data.
ampl();
996 if (!info_list.empty())
997 hitInfo = std::max_element(info_list.begin(), info_list.end())->
second.get();
1002 info.ot_bit = signalInElectrons > theHIPThresholdInE ?
true :
false;
1005 float charge_frac =
l.first / signalInElectrons;
1007 info.simInfoList.push_back({charge_frac,
l.second.get()});
1010 digi_map.insert({
s.first,
info});
1020 const int max_limit = 10;
1032 if (temp_signal > kink_point)
1033 temp_signal = std::floor((temp_signal - kink_point) / (
pow(2, dualslope_param - 1))) + kink_point;
1036 LogTrace(
"Phase2TrackerDigitizerAlgorithm")
1037 <<
" DetId " << detID <<
" signal_in_elec " << signal_in_elec <<
" threshold " <<
threshold
1038 <<
" signal_above_thr " << signal_in_elec -
threshold <<
" temp conversion "
1040 << temp_signal <<
" signal_in_adc " << signal_in_adc;
1042 return signal_in_adc;