Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-04-18 22:24:06

0001 #include <typeinfo>
0002 #include <iostream>
0003 #include <cmath>
0004 
0005 #include "SimTracker/SiPhase2Digitizer/plugins/Phase2TrackerDigitizerAlgorithm.h"
0006 
0007 #include "SimGeneral/NoiseGenerators/interface/GaussianTailNoiseGenerator.h"
0008 #include "SimTracker/Common/interface/SiG4UniversalFluctuation.h"
0009 
0010 #include "CLHEP/Random/RandGaussQ.h"
0011 
0012 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0013 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0014 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0015 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0016 
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019 #include "FWCore/ServiceRegistry/interface/Service.h"
0020 #include "FWCore/Utilities/interface/Exception.h"
0021 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationOfflineSimService.h"
0022 
0023 #include "DataFormats/DetId/interface/DetId.h"
0024 #include "CondFormats/SiPixelObjects/interface/GlobalPixel.h"
0025 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
0026 #include "CondFormats/SiPixelObjects/interface/LocalPixel.h"
0027 #include "CondFormats/SiPhase2TrackerObjects/interface/SiPhase2OuterTrackerLorentzAngle.h"
0028 
0029 // Geometry
0030 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0031 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0032 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0033 
0034 using namespace edm;
0035 
0036 namespace ph2tkdigialgo {
0037   // Mass in MeV
0038   constexpr double m_pion = 139.571;
0039   constexpr double m_kaon = 493.677;
0040   constexpr double m_electron = 0.511;
0041   constexpr double m_muon = 105.658;
0042   constexpr double m_proton = 938.272;
0043 }  // namespace ph2tkdigialgo
0044 Phase2TrackerDigitizerAlgorithm::Phase2TrackerDigitizerAlgorithm(const edm::ParameterSet& conf_common,
0045                                                                  const edm::ParameterSet& conf_specific,
0046                                                                  edm::ConsumesCollector iC)
0047     : _signal(),
0048       makeDigiSimLinks_(conf_common.getUntrackedParameter<bool>("makeDigiSimLinks", true)),
0049       use_ineff_from_db_(conf_specific.getParameter<bool>("Inefficiency_DB")),
0050       use_module_killing_(conf_specific.getParameter<bool>("KillModules")),    // boolean to kill or not modules
0051       use_deadmodule_DB_(conf_specific.getParameter<bool>("DeadModules_DB")),  // boolean to access dead modules from DB
0052       // boolean to access Lorentz angle from DB
0053       use_LorentzAngle_DB_(conf_specific.getParameter<bool>("LorentzAngle_DB")),
0054 
0055       // get dead module from cfg file
0056       deadModules_(use_deadmodule_DB_ ? Parameters() : conf_specific.getParameter<Parameters>("DeadModules")),
0057 
0058       // Common pixel parameters
0059       // These are parameters which are not likely to be changed
0060       GeVperElectron_(3.61E-09),                                      // 1 electron(3.61eV, 1keV(277e, mod 9/06 d.k.
0061       alpha2Order_(conf_specific.getParameter<bool>("Alpha2Order")),  // switch on/off of E.B effect
0062       addXtalk_(conf_specific.getParameter<bool>("AddXTalk")),
0063       // Interstrip Coupling - Not used in PixelDigitizerAlgorithm
0064       interstripCoupling_(conf_specific.getParameter<double>("InterstripCoupling")),
0065 
0066       Sigma0_(conf_specific.getParameter<double>("SigmaZero")),       // Charge diffusion constant 7->3.7
0067       SigmaCoeff_(conf_specific.getParameter<double>("SigmaCoeff")),  // delta in the diffusion across the strip pitch
0068       // (set between 0 to 0.9,  0-->flat Sigma0, 1-->Sigma_min=0 & Sigma_max=2*Sigma0
0069       // D.B.: Dist300 replaced by moduleThickness, may not work with partially depleted sensors but works otherwise
0070       // Dist300(0.0300),                                          //   normalized to 300micron Silicon
0071 
0072       // Charge integration spread on the collection plane
0073       clusterWidth_(conf_specific.getParameter<double>("ClusterWidth")),
0074 
0075       // Allowed modes of readout which has following values :
0076       // 0          ---> Digital or binary readout
0077       // -1         ---> Analog readout, current digitizer (Inner Pixel) (TDR version) with no threshold subtraction
0078       // Analog readout with dual slope with the "second" slope being 1/2^(n-1) and threshold subtraction (n = 1, 2, 3,4)
0079       thePhase2ReadoutMode_(conf_specific.getParameter<int>("Phase2ReadoutMode")),
0080 
0081       // ADC calibration 1adc count(135e.
0082       // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
0083       // Be careful, this parameter is also used in SiPixelDet.cc to
0084       // calculate the noise in adc counts from noise in electrons.
0085       // Both defaults should be the same.
0086       theElectronPerADC_(conf_specific.getParameter<double>("ElectronPerAdc")),
0087 
0088       // ADC saturation value, 255(8bit adc.
0089       theAdcFullScale_(conf_specific.getParameter<int>("AdcFullScale")),
0090 
0091       // Noise in electrons:
0092       // Pixel cell noise, relevant for generating noisy pixels
0093       theNoiseInElectrons_(conf_specific.getParameter<double>("NoiseInElectrons")),
0094 
0095       // Fill readout noise, including all readout chain, relevant for smearing
0096       theReadoutNoise_(conf_specific.getParameter<double>("ReadoutNoiseInElec")),
0097 
0098       // Threshold in units of noise:
0099       // thePixelThreshold(conf.getParameter<double>("ThresholdInNoiseUnits")),
0100       // Pixel threshold in electron units.
0101       theThresholdInE_Endcap_(conf_specific.getParameter<double>("ThresholdInElectrons_Endcap")),
0102       theThresholdInE_Barrel_(conf_specific.getParameter<double>("ThresholdInElectrons_Barrel")),
0103 
0104       // Add threshold gaussian smearing:
0105       theThresholdSmearing_Endcap_(conf_specific.getParameter<double>("ThresholdSmearing_Endcap")),
0106       theThresholdSmearing_Barrel_(conf_specific.getParameter<double>("ThresholdSmearing_Barrel")),
0107 
0108       // Add HIP Threshold in electron units.
0109       theHIPThresholdInE_Endcap_(conf_specific.getParameter<double>("HIPThresholdInElectrons_Endcap")),
0110       theHIPThresholdInE_Barrel_(conf_specific.getParameter<double>("HIPThresholdInElectrons_Barrel")),
0111 
0112       // theTofCut 12.5, cut in particle TOD +/- 12.5ns
0113       theTofLowerCut_(conf_specific.getParameter<double>("TofLowerCut")),
0114       theTofUpperCut_(conf_specific.getParameter<double>("TofUpperCut")),
0115 
0116       // Get the Lorentz angle from the cfg file:
0117       tanLorentzAnglePerTesla_Endcap_(
0118           use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Endcap")),
0119       tanLorentzAnglePerTesla_Barrel_(
0120           use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Barrel")),
0121 
0122       // Add noise
0123       addNoise_(conf_specific.getParameter<bool>("AddNoise")),
0124 
0125       // Add noisy pixels
0126       addNoisyPixels_(conf_specific.getParameter<bool>("AddNoisyPixels")),
0127 
0128       // Fluctuate charge in track subsegments
0129       fluctuateCharge_(conf_specific.getUntrackedParameter<bool>("FluctuateCharge", true)),
0130 
0131       // Control the pixel inefficiency
0132       addPixelInefficiency_(conf_specific.getParameter<bool>("AddInefficiency")),
0133 
0134       // Add threshold gaussian smearing:
0135       addThresholdSmearing_(conf_specific.getParameter<bool>("AddThresholdSmearing")),
0136 
0137       // Add some pseudo-red damage
0138       pseudoRadDamage_(conf_specific.exists("PseudoRadDamage") ? conf_specific.getParameter<double>("PseudoRadDamage")
0139                                                                : double(0.0)),
0140       pseudoRadDamageRadius_(conf_specific.exists("PseudoRadDamageRadius")
0141                                  ? conf_specific.getParameter<double>("PseudoRadDamageRadius")
0142                                  : double(0.0)),
0143 
0144       // delta cutoff in MeV, has to be same as in OSCAR(0.030/cmsim=1.0 MeV
0145       // tMax(0.030), // In MeV.
0146       // tMax(conf.getUntrackedParameter<double>("DeltaProductionCut",0.030)),
0147       tMax_(conf_common.getParameter<double>("DeltaProductionCut")),
0148 
0149       badPixels_(conf_specific.getParameter<Parameters>("CellsToKill")),
0150 
0151       fluctuate_(fluctuateCharge_ ? std::make_unique<SiG4UniversalFluctuation>() : nullptr),
0152       theNoiser_(addNoise_ ? std::make_unique<GaussianTailNoiseGenerator>() : nullptr),
0153       theSiPixelGainCalibrationService_(
0154           use_ineff_from_db_ ? std::make_unique<SiPixelGainCalibrationOfflineSimService>(conf_specific, iC) : nullptr),
0155       subdetEfficiencies_(conf_specific) {
0156   LogDebug("Phase2TrackerDigitizerAlgorithm")
0157       << "Phase2TrackerDigitizerAlgorithm constructed\n"
0158       << "Configuration parameters:\n"
0159       << "Threshold/Gain = "
0160       << "threshold in electron Endcap = " << theThresholdInE_Endcap_
0161       << "\nthreshold in electron Barrel = " << theThresholdInE_Barrel_ << " ElectronPerADC " << theElectronPerADC_
0162       << " ADC Scale (in bits) " << theAdcFullScale_ << " The delta cut-off is set to " << tMax_ << " pix-inefficiency "
0163       << addPixelInefficiency_;
0164 }
0165 
0166 Phase2TrackerDigitizerAlgorithm::~Phase2TrackerDigitizerAlgorithm() {
0167   LogDebug("Phase2TrackerDigitizerAlgorithm") << "Phase2TrackerDigitizerAlgorithm deleted";
0168 }
0169 
0170 Phase2TrackerDigitizerAlgorithm::SubdetEfficiencies::SubdetEfficiencies(const edm::ParameterSet& conf) {
0171   barrel_efficiencies = conf.getParameter<std::vector<double> >("EfficiencyFactors_Barrel");
0172   endcap_efficiencies = conf.getParameter<std::vector<double> >("EfficiencyFactors_Endcap");
0173 }
0174 void Phase2TrackerDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
0175                                                         std::vector<PSimHit>::const_iterator inputEnd,
0176                                                         const size_t inputBeginGlobalIndex,
0177                                                         const uint32_t tofBin,
0178                                                         const Phase2TrackerGeomDetUnit* pixdet,
0179                                                         const GlobalVector& bfield) {
0180   // produce SignalPoint's for all SimHit's in detector
0181   // Loop over hits
0182   uint32_t detId = pixdet->geographicalId().rawId();
0183   size_t simHitGlobalIndex = inputBeginGlobalIndex;  // This needs to be stored to create the digi-sim link later
0184 
0185   // find the relevant hits
0186   std::vector<PSimHit> matchedSimHits;
0187   std::copy_if(inputBegin, inputEnd, std::back_inserter(matchedSimHits), [detId](auto const& hit) -> bool {
0188     return hit.detUnitId() == detId;
0189   });
0190   // loop over a much reduced set of SimHits
0191   for (auto const& hit : matchedSimHits) {
0192     LogDebug("Phase2DigitizerAlgorithm") << hit.particleType() << " " << hit.pabs() << " " << hit.energyLoss() << " "
0193                                          << hit.tof() << " " << hit.trackId() << " " << hit.processType() << " "
0194                                          << hit.detUnitId() << hit.entryPoint() << " " << hit.exitPoint();
0195     double signalScale = 1.0;
0196     // fill collection_points for this SimHit, indpendent of topology
0197     if (select_hit(hit, (pixdet->surface().toGlobal(hit.localPosition()).mag() * c_inv), signalScale)) {
0198       const auto& ionization_points = primary_ionization(hit);  // fills ionization_points
0199 
0200       // transforms ionization_points -> collection_points
0201       const auto& collection_points = drift(hit, pixdet, bfield, ionization_points);
0202 
0203       // compute induced signal on readout elements and add to _signal
0204       // hit needed only for SimHit<-->Digi link
0205       induce_signal(hit, simHitGlobalIndex, tofBin, pixdet, collection_points);
0206     }
0207     ++simHitGlobalIndex;
0208   }
0209 }
0210 // =================================================================
0211 //
0212 // Generate primary ionization along the track segment.
0213 // Divide the track into small sub-segments
0214 //
0215 // =================================================================
0216 std::vector<DigitizerUtility::EnergyDepositUnit> Phase2TrackerDigitizerAlgorithm::primary_ionization(
0217     const PSimHit& hit) const {
0218   // Straight line approximation for trajectory inside active media
0219   constexpr float SegmentLength = 0.0010;  // in cm (10 microns)
0220   // Get the 3D segment direction vector
0221   LocalVector direction = hit.exitPoint() - hit.entryPoint();
0222 
0223   float eLoss = hit.energyLoss();  // Eloss in GeV
0224   float length = direction.mag();  // Track length in Silicon
0225 
0226   int NumberOfSegments = static_cast<int>(length / SegmentLength);  // Number of segments
0227   if (NumberOfSegments < 1)
0228     NumberOfSegments = 1;
0229   LogDebug("Phase2TrackerDigitizerAlgorithm")
0230       << "enter primary_ionzation " << NumberOfSegments << " shift = " << hit.exitPoint().x() - hit.entryPoint().x()
0231       << " " << hit.exitPoint().y() - hit.entryPoint().y() << " " << hit.exitPoint().z() - hit.entryPoint().z() << " "
0232       << hit.particleType() << " " << hit.pabs();
0233 
0234   std::vector<float> elossVector;
0235   elossVector.reserve(NumberOfSegments);
0236   if (fluctuateCharge_) {
0237     // Generate fluctuated charge points
0238     elossVector = fluctuateEloss(hit.particleType(), hit.pabs(), eLoss, length, NumberOfSegments);
0239   } else {
0240     float averageEloss = eLoss / NumberOfSegments;
0241     elossVector.resize(NumberOfSegments, averageEloss);
0242   }
0243 
0244   std::vector<DigitizerUtility::EnergyDepositUnit> ionization_points;
0245   ionization_points.reserve(NumberOfSegments);  // set size
0246   // loop over segments
0247   for (size_t i = 0; i < elossVector.size(); ++i) {
0248     // Divide the segment into equal length subsegments
0249     Local3DPoint point = hit.entryPoint() + ((i + 0.5) / NumberOfSegments) * direction;
0250     float energy = elossVector[i] / GeVperElectron_;  // Convert charge to elec.
0251 
0252     DigitizerUtility::EnergyDepositUnit edu(energy, point);  // define position,energy point
0253     ionization_points.push_back(edu);                        // save
0254     LogDebug("Phase2TrackerDigitizerAlgorithm")
0255         << i << " " << edu.x() << " " << edu.y() << " " << edu.z() << " " << edu.energy();
0256   }
0257   return ionization_points;
0258 }
0259 //==============================================================================
0260 //
0261 // Fluctuate the charge comming from a small (10um) track segment.
0262 // Use the G4 routine. For mip pions for the moment.
0263 //
0264 //==============================================================================
0265 std::vector<float> Phase2TrackerDigitizerAlgorithm::fluctuateEloss(
0266     int pid, float particleMomentum, float eloss, float length, int NumberOfSegs) const {
0267   double particleMass = ph2tkdigialgo::m_pion;  // Mass in MeV, assume pion
0268   pid = std::abs(pid);
0269   if (pid != 211) {  // Mass in MeV
0270     if (pid == 11)
0271       particleMass = ph2tkdigialgo::m_electron;
0272     else if (pid == 13)
0273       particleMass = ph2tkdigialgo::m_muon;
0274     else if (pid == 321)
0275       particleMass = ph2tkdigialgo::m_kaon;
0276     else if (pid == 2212)
0277       particleMass = ph2tkdigialgo::m_proton;
0278   }
0279   // What is the track segment length.
0280   float segmentLength = length / NumberOfSegs;
0281 
0282   // Generate charge fluctuations.
0283   float sum = 0.;
0284   double segmentEloss = (1000. * eloss) / NumberOfSegs;  //eloss in MeV
0285   std::vector<float> elossVector;
0286   elossVector.reserve(NumberOfSegs);
0287   for (int i = 0; i < NumberOfSegs; ++i) {
0288     // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
0289     // track segment length in mm, segment eloss in MeV
0290     // Returns fluctuated eloss in MeV
0291     double deltaCutoff = tMax_;  // the cutoff is sometimes redefined inside, so fix it.
0292     float de = fluctuate_->SampleFluctuations(particleMomentum * 1000.,
0293                                               particleMass,
0294                                               deltaCutoff,
0295                                               segmentLength * 10.,
0296                                               segmentEloss,
0297                                               rengine_) /
0298                1000.;  //convert to GeV
0299     elossVector.push_back(de);
0300     sum += de;
0301   }
0302   if (sum > 0.) {  // if fluctuations give eloss>0.
0303     // Rescale to the same total eloss
0304     float ratio = eloss / sum;
0305     std::transform(
0306         std::begin(elossVector), std::end(elossVector), std::begin(elossVector), [&ratio](auto const& c) -> float {
0307           return c * ratio;
0308         });  // use a simple lambda expression
0309   } else {   // if fluctuations gives 0 eloss
0310     float averageEloss = eloss / NumberOfSegs;
0311     elossVector.resize(NumberOfSegs, averageEloss);
0312   }
0313   return elossVector;
0314 }
0315 
0316 // ======================================================================
0317 //
0318 // Drift the charge segments to the sensor surface (collection plane)
0319 // Include the effect of E-field and B-field
0320 //
0321 // =====================================================================
0322 std::vector<DigitizerUtility::SignalPoint> Phase2TrackerDigitizerAlgorithm::drift(
0323     const PSimHit& hit,
0324     const Phase2TrackerGeomDetUnit* pixdet,
0325     const GlobalVector& bfield,
0326     const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points) const {
0327   LogDebug("Phase2TrackerDigitizerAlgorithm") << "enter drift ";
0328 
0329   std::vector<DigitizerUtility::SignalPoint> collection_points;
0330   collection_points.reserve(ionization_points.size());                     // set size
0331   LocalVector driftDir = DriftDirection(pixdet, bfield, hit.detUnitId());  // get the charge drift direction
0332   if (driftDir.z() == 0.) {
0333     LogWarning("Phase2TrackerDigitizerAlgorithm") << " pxlx: drift in z is zero ";
0334     return collection_points;
0335   }
0336 
0337   float TanLorenzAngleX = driftDir.x();                                       // tangent of Lorentz angle
0338   float TanLorenzAngleY = 0.;                                                 // force to 0, driftDir.y()/driftDir.z();
0339   float dir_z = driftDir.z();                                                 // The z drift direction
0340   float CosLorenzAngleX = 1. / std::sqrt(1. + std::pow(TanLorenzAngleX, 2));  // cosine to estimate the path length
0341   float CosLorenzAngleY = 1.;
0342   if (alpha2Order_) {
0343     TanLorenzAngleY = driftDir.y();
0344     CosLorenzAngleY = 1. / std::sqrt(1. + std::pow(TanLorenzAngleY, 2));  // cosine
0345   }
0346 
0347   float moduleThickness = pixdet->specificSurface().bounds().thickness();
0348   float stripPitch = pixdet->specificTopology().pitch().first;
0349 
0350   LogDebug("Phase2TrackerDigitizerAlgorithm")
0351       << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY << " " << CosLorenzAngleX << " "
0352       << CosLorenzAngleY << " " << moduleThickness * TanLorenzAngleX << " " << driftDir;
0353 
0354   for (auto const& val : ionization_points) {
0355     // position
0356     float SegX = val.x(), SegY = val.y(), SegZ = val.z();
0357 
0358     // Distance from the collection plane
0359     // DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
0360     // Include explixitely the E drift direction (for CMS dir_z=-1)
0361 
0362     // Distance between charge generation and collection
0363     float driftDistance = moduleThickness / 2. - (dir_z * SegZ);  // Drift to -z
0364 
0365     if (driftDistance < 0.)
0366       driftDistance = 0.;
0367     else if (driftDistance > moduleThickness)
0368       driftDistance = moduleThickness;
0369 
0370     // Assume full depletion now, partial depletion will come later.
0371     float XDriftDueToMagField = driftDistance * TanLorenzAngleX;
0372     float YDriftDueToMagField = driftDistance * TanLorenzAngleY;
0373 
0374     // Shift cloud center
0375     float CloudCenterX = SegX + XDriftDueToMagField;
0376     float CloudCenterY = SegY + YDriftDueToMagField;
0377 
0378     // Calculate how long is the charge drift path
0379     // Actual Drift Lentgh
0380     float driftLength =
0381         std::sqrt(std::pow(driftDistance, 2) + std::pow(XDriftDueToMagField, 2) + std::pow(YDriftDueToMagField, 2));
0382 
0383     // What is the charge diffusion after this path
0384     // Sigma0=0.00037 is for 300um thickness (make sure moduleThickness is in [cm])
0385     float Sigma = std::sqrt(driftLength / moduleThickness) * Sigma0_ * moduleThickness / 0.0300;
0386     // D.B.: sigmaCoeff=0 means no modulation
0387     if (SigmaCoeff_)
0388       Sigma *= (SigmaCoeff_ * std::pow(cos(SegX * M_PI / stripPitch), 2) + 1);
0389     // NB: divided by 4 to get a periodicity of stripPitch
0390 
0391     // Project the diffusion sigma on the collection plane
0392     float Sigma_x = Sigma / CosLorenzAngleX;
0393     float Sigma_y = Sigma / CosLorenzAngleY;
0394 
0395     // Insert a charge loss due to Rad Damage here
0396     float energyOnCollector = val.energy();  // The energy that reaches the collector
0397 
0398     // pseudoRadDamage
0399     if (pseudoRadDamage_) {
0400       float moduleRadius = pixdet->surface().position().perp();
0401       if (moduleRadius <= pseudoRadDamageRadius_) {
0402         float kValue = pseudoRadDamage_ / std::pow(moduleRadius, 2);
0403         energyOnCollector *= exp(-1 * kValue * driftDistance / moduleThickness);
0404       }
0405     }
0406     LogDebug("Phase2TrackerDigitizerAlgorithm")
0407         << "Dift DistanceZ = " << driftDistance << " module thickness = " << moduleThickness
0408         << " Start Energy = " << val.energy() << " Energy after loss= " << energyOnCollector;
0409     DigitizerUtility::SignalPoint sp(CloudCenterX, CloudCenterY, Sigma_x, Sigma_y, hit.tof(), energyOnCollector);
0410 
0411     // Load the Charge distribution parameters
0412     collection_points.push_back(sp);
0413   }
0414   return collection_points;
0415 }
0416 
0417 // ====================================================================
0418 //
0419 // Induce the signal on the collection plane of the active sensor area.
0420 void Phase2TrackerDigitizerAlgorithm::induce_signal(
0421     const PSimHit& hit,
0422     const size_t hitIndex,
0423     const uint32_t tofBin,
0424     const Phase2TrackerGeomDetUnit* pixdet,
0425     const std::vector<DigitizerUtility::SignalPoint>& collection_points) {
0426   // X  - Rows, Left-Right, 160, (1.6cm)   for barrel
0427   // Y  - Columns, Down-Up, 416, (6.4cm)
0428   const Phase2TrackerTopology* topol = &pixdet->specificTopology();
0429   uint32_t detID = pixdet->geographicalId().rawId();
0430   signal_map_type& theSignal = _signal[detID];
0431 
0432   LogDebug("Phase2TrackerDigitizerAlgorithm")
0433       << " enter induce_signal, " << topol->pitch().first << " " << topol->pitch().second;
0434 
0435   // local map to store pixels hit by 1 Hit.
0436   using hit_map_type = std::map<int, float, std::less<int> >;
0437   hit_map_type hit_signal;
0438 
0439   // Assign signals to readout channels and store sorted by channel number
0440   // Iterate over collection points on the collection plane
0441   for (auto const& v : collection_points) {
0442     float CloudCenterX = v.position().x();  // Charge position in x
0443     float CloudCenterY = v.position().y();  //                 in y
0444     float SigmaX = v.sigma_x();             // Charge spread in x
0445     float SigmaY = v.sigma_y();             //               in y
0446     float Charge = v.amplitude();           // Charge amplitude
0447 
0448     LogDebug("Phase2TrackerDigitizerAlgorithm") << " cloud " << v.position().x() << " " << v.position().y() << " "
0449                                                 << v.sigma_x() << " " << v.sigma_y() << " " << v.amplitude();
0450 
0451     // Find the maximum cloud spread in 2D plane , assume 3*sigma
0452     float CloudRight = CloudCenterX + clusterWidth_ * SigmaX;
0453     float CloudLeft = CloudCenterX - clusterWidth_ * SigmaX;
0454     float CloudUp = CloudCenterY + clusterWidth_ * SigmaY;
0455     float CloudDown = CloudCenterY - clusterWidth_ * SigmaY;
0456 
0457     // Define 2D cloud limit points
0458     LocalPoint PointRightUp = LocalPoint(CloudRight, CloudUp);
0459     LocalPoint PointLeftDown = LocalPoint(CloudLeft, CloudDown);
0460 
0461     // This points can be located outside the sensor area.
0462     // The conversion to measurement point does not check for that
0463     // so the returned pixel index might be wrong (outside range).
0464     // We rely on the limits check below to fix this.
0465     // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
0466     // PIXEL ARE IS LOST, it should not be collected.
0467 
0468     // Convert the 2D points to pixel indices
0469     MeasurementPoint mp = topol->measurementPosition(PointRightUp);
0470     int IPixRightUpX = static_cast<int>(std::floor(mp.x()));  // cast reqd.
0471     int IPixRightUpY = static_cast<int>(std::floor(mp.y()));
0472     LogDebug("Phase2TrackerDigitizerAlgorithm")
0473         << " right-up " << PointRightUp << " " << mp.x() << " " << mp.y() << " " << IPixRightUpX << " " << IPixRightUpY;
0474 
0475     mp = topol->measurementPosition(PointLeftDown);
0476     int IPixLeftDownX = static_cast<int>(std::floor(mp.x()));
0477     int IPixLeftDownY = static_cast<int>(std::floor(mp.y()));
0478     LogDebug("Phase2TrackerDigitizerAlgorithm") << " left-down " << PointLeftDown << " " << mp.x() << " " << mp.y()
0479                                                 << " " << IPixLeftDownX << " " << IPixLeftDownY;
0480 
0481     // Check detector limits to correct for pixels outside range.
0482     int numColumns = topol->ncolumns();  // det module number of cols&rows
0483     int numRows = topol->nrows();
0484 
0485     IPixRightUpX = numRows > IPixRightUpX ? IPixRightUpX : numRows - 1;
0486     IPixRightUpY = numColumns > IPixRightUpY ? IPixRightUpY : numColumns - 1;
0487     IPixLeftDownX = 0 < IPixLeftDownX ? IPixLeftDownX : 0;
0488     IPixLeftDownY = 0 < IPixLeftDownY ? IPixLeftDownY : 0;
0489 
0490     // First integrate charge strips in x
0491     hit_map_type x;
0492     for (int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) {  // loop over x index
0493       float xLB, LowerBound;
0494       // Why is set to 0 if ix=0, does it meen that we accept charge
0495       // outside the sensor?
0496       if (ix == 0 || SigmaX == 0.) {  // skip for surface segemnts
0497         LowerBound = 0.;
0498       } else {
0499         mp = MeasurementPoint(ix, 0.0);
0500         xLB = topol->localPosition(mp).x();
0501         LowerBound = 1 - calcQ((xLB - CloudCenterX) / SigmaX);
0502       }
0503 
0504       float xUB, UpperBound;
0505       if (ix == numRows - 1 || SigmaX == 0.) {
0506         UpperBound = 1.;
0507       } else {
0508         mp = MeasurementPoint(ix + 1, 0.0);
0509         xUB = topol->localPosition(mp).x();
0510         UpperBound = 1. - calcQ((xUB - CloudCenterX) / SigmaX);
0511       }
0512       float TotalIntegrationRange = UpperBound - LowerBound;  // get strip
0513       x.emplace(ix, TotalIntegrationRange);                   // save strip integral
0514     }
0515 
0516     // Now integrate strips in y
0517     hit_map_type y;
0518     for (int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) {  // loop over y index
0519       float yLB, LowerBound;
0520       if (iy == 0 || SigmaY == 0.) {
0521         LowerBound = 0.;
0522       } else {
0523         mp = MeasurementPoint(0.0, iy);
0524         yLB = topol->localPosition(mp).y();
0525         LowerBound = 1. - calcQ((yLB - CloudCenterY) / SigmaY);
0526       }
0527 
0528       float yUB, UpperBound;
0529       if (iy == numColumns - 1 || SigmaY == 0.) {
0530         UpperBound = 1.;
0531       } else {
0532         mp = MeasurementPoint(0.0, iy + 1);
0533         yUB = topol->localPosition(mp).y();
0534         UpperBound = 1. - calcQ((yUB - CloudCenterY) / SigmaY);
0535       }
0536 
0537       float TotalIntegrationRange = UpperBound - LowerBound;
0538       y.emplace(iy, TotalIntegrationRange);  // save strip integral
0539     }
0540 
0541     // Get the 2D charge integrals by folding x and y strips
0542     for (int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) {    // loop over x index
0543       for (int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) {  // loop over y index
0544         float ChargeFraction = Charge * x[ix] * y[iy];
0545         int chanFired = -1;
0546         if (ChargeFraction > 0.) {
0547           chanFired =
0548               pixelFlag_ ? PixelDigi::pixelToChannel(ix, iy) : Phase2TrackerDigi::pixelToChannel(ix, iy);  // Get index
0549           // Load the amplitude
0550           hit_signal[chanFired] += ChargeFraction;
0551         }
0552 
0553         mp = MeasurementPoint(ix, iy);
0554         LocalPoint lp = topol->localPosition(mp);
0555         int chan = topol->channel(lp);
0556 
0557         LogDebug("Phase2TrackerDigitizerAlgorithm")
0558             << " pixel " << ix << " " << iy << " - "
0559             << " " << chanFired << " " << ChargeFraction << " " << mp.x() << " " << mp.y() << " " << lp.x() << " "
0560             << lp.y() << " "  // givex edge position
0561             << chan;          // edge belongs to previous ?
0562       }
0563     }
0564   }
0565   // Fill the global map with all hit pixels from this event
0566   float corr_time = hit.tof() - pixdet->surface().toGlobal(hit.localPosition()).mag() * c_inv;
0567   for (auto const& hit_s : hit_signal) {
0568     int chan = hit_s.first;
0569     theSignal[chan] +=
0570         (makeDigiSimLinks_ ? DigitizerUtility::Amplitude(hit_s.second, &hit, hit_s.second, corr_time, hitIndex, tofBin)
0571                            : DigitizerUtility::Amplitude(hit_s.second, nullptr, hit_s.second));
0572   }
0573 }
0574 // ======================================================================
0575 //
0576 //  Add electronic noise to pixel charge
0577 //
0578 // ======================================================================
0579 void Phase2TrackerDigitizerAlgorithm::add_noise(const Phase2TrackerGeomDetUnit* pixdet) {
0580   uint32_t detID = pixdet->geographicalId().rawId();
0581   signal_map_type& theSignal = _signal[detID];
0582   for (auto& s : theSignal) {
0583     float noise = gaussDistribution_->fire();
0584     if ((s.second.ampl() + noise) < 0.)
0585       s.second.set(0);
0586     else
0587       s.second += noise;
0588   }
0589 }
0590 // ======================================================================
0591 //
0592 //  Add  Cross-talk contribution
0593 //
0594 // ======================================================================
0595 void Phase2TrackerDigitizerAlgorithm::add_cross_talk(const Phase2TrackerGeomDetUnit* pixdet) {
0596   uint32_t detID = pixdet->geographicalId().rawId();
0597   signal_map_type& theSignal = _signal[detID];
0598   signal_map_type signalNew;
0599   const Phase2TrackerTopology* topol = &pixdet->specificTopology();
0600   int numRows = topol->nrows();
0601 
0602   for (auto& s : theSignal) {
0603     float signalInElectrons = s.second.ampl();  // signal in electrons
0604 
0605     std::pair<int, int> hitChan;
0606     if (pixelFlag_)
0607       hitChan = PixelDigi::channelToPixel(s.first);
0608     else
0609       hitChan = Phase2TrackerDigi::channelToPixel(s.first);
0610 
0611     float signalInElectrons_Xtalk = signalInElectrons * interstripCoupling_;
0612     // subtract the charge which will be shared
0613     s.second.set(signalInElectrons - signalInElectrons_Xtalk);
0614 
0615     if (hitChan.first != 0) {
0616       auto XtalkPrev = std::make_pair(hitChan.first - 1, hitChan.second);
0617       int chanXtalkPrev = pixelFlag_ ? PixelDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second)
0618                                      : Phase2TrackerDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second);
0619       signalNew.emplace(chanXtalkPrev, DigitizerUtility::Amplitude(signalInElectrons_Xtalk, nullptr, -1.0));
0620     }
0621     if (hitChan.first < numRows - 1) {
0622       auto XtalkNext = std::make_pair(hitChan.first + 1, hitChan.second);
0623       int chanXtalkNext = pixelFlag_ ? PixelDigi::pixelToChannel(XtalkNext.first, XtalkNext.second)
0624                                      : Phase2TrackerDigi::pixelToChannel(XtalkNext.first, XtalkNext.second);
0625       signalNew.emplace(chanXtalkNext, DigitizerUtility::Amplitude(signalInElectrons_Xtalk, nullptr, -1.0));
0626     }
0627   }
0628   for (auto const& l : signalNew) {
0629     int chan = l.first;
0630     auto iter = theSignal.find(chan);
0631     if (iter != theSignal.end()) {
0632       theSignal[chan] += l.second.ampl();
0633     } else {
0634       theSignal.emplace(chan, DigitizerUtility::Amplitude(l.second.ampl(), nullptr, -1.0));
0635     }
0636   }
0637 }
0638 
0639 // ======================================================================
0640 //
0641 //  Add noise on non-hit cells
0642 //
0643 // ======================================================================
0644 void Phase2TrackerDigitizerAlgorithm::add_noisy_cells(const Phase2TrackerGeomDetUnit* pixdet, float thePixelThreshold) {
0645   uint32_t detID = pixdet->geographicalId().rawId();
0646   signal_map_type& theSignal = _signal[detID];
0647   const Phase2TrackerTopology* topol = &pixdet->specificTopology();
0648 
0649   int numColumns = topol->ncolumns();  // det module number of cols&rows
0650   int numRows = topol->nrows();
0651 
0652   int numberOfPixels = numRows * numColumns;
0653   std::map<int, float, std::less<int> > otherPixels;
0654 
0655   theNoiser_->generate(numberOfPixels,
0656                        thePixelThreshold,     //thr. in un. of nois
0657                        theNoiseInElectrons_,  // noise in elec.
0658                        otherPixels,
0659                        rengine_);
0660 
0661   LogDebug("Phase2TrackerDigitizerAlgorithm")
0662       << " Add noisy pixels " << numRows << " " << numColumns << " " << theNoiseInElectrons_ << " "
0663       << theThresholdInE_Endcap_ << "  " << theThresholdInE_Barrel_ << " " << numberOfPixels << " "
0664       << otherPixels.size();
0665 
0666   // Add noisy pixels
0667   for (auto const& el : otherPixels) {
0668     int iy = el.first / numRows;
0669     if (iy < 0 || iy > numColumns - 1)
0670       LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in iy " << iy;
0671 
0672     int ix = el.first - iy * numRows;
0673     if (ix < 0 || ix > numRows - 1)
0674       LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in ix " << ix;
0675 
0676     int chan = pixelFlag_ ? PixelDigi::pixelToChannel(ix, iy) : Phase2TrackerDigi::pixelToChannel(ix, iy);
0677 
0678     LogDebug("Phase2TrackerDigitizerAlgorithm")
0679         << " Storing noise = " << el.first << " " << el.second << " " << ix << " " << iy << " " << chan;
0680 
0681     if (theSignal[chan] == 0)
0682       theSignal[chan] = DigitizerUtility::Amplitude(el.second, nullptr, -1.);
0683   }
0684 }
0685 // ============================================================================
0686 //
0687 // Simulate the readout inefficiencies.
0688 // Delete a selected number of single pixels, dcols and rocs.
0689 void Phase2TrackerDigitizerAlgorithm::pixel_inefficiency(const SubdetEfficiencies& eff,
0690                                                          const Phase2TrackerGeomDetUnit* pixdet,
0691                                                          const TrackerTopology* tTopo) {
0692   uint32_t detID = pixdet->geographicalId().rawId();
0693   signal_map_type& theSignal = _signal[detID];  // check validity
0694 
0695   // Predefined efficiencies
0696   float subdetEfficiency = 1.0;
0697 
0698   // setup the chip indices conversion
0699   uint32_t Subid = DetId(detID).subdetId();
0700   if (Subid == PixelSubdetector::PixelBarrel || Subid == StripSubdetector::TOB) {  // barrel layers
0701     uint32_t layerIndex = tTopo->pxbLayer(detID);
0702     if (layerIndex - 1 < eff.barrel_efficiencies.size())
0703       subdetEfficiency = eff.barrel_efficiencies[layerIndex - 1];
0704   } else {  // forward disks
0705     uint32_t diskIndex = 2 * tTopo->pxfDisk(detID) - tTopo->pxfSide(detID);
0706     if (diskIndex - 1 < eff.endcap_efficiencies.size())
0707       subdetEfficiency = eff.endcap_efficiencies[diskIndex - 1];
0708   }
0709 
0710   LogDebug("Phase2TrackerDigitizerAlgorithm") << " enter pixel_inefficiency " << subdetEfficiency;
0711 
0712   // Now loop again over pixels to kill some of them.
0713   // Loop over hits, amplitude in electrons, channel = coded row,col
0714   for (auto& s : theSignal) {
0715     float rand = rengine_->flat();
0716     if (rand > subdetEfficiency) {
0717       // make amplitude =0
0718       s.second.set(0.);  // reset amplitude
0719     }
0720   }
0721 }
0722 void Phase2TrackerDigitizerAlgorithm::initializeEvent(CLHEP::HepRandomEngine& eng) {
0723   if (addNoise_ || addPixelInefficiency_ || fluctuateCharge_ || addThresholdSmearing_) {
0724     gaussDistribution_ = std::make_unique<CLHEP::RandGaussQ>(eng, 0., theReadoutNoise_);
0725   }
0726   // Threshold smearing with gaussian distribution:
0727   if (addThresholdSmearing_) {
0728     smearedThreshold_Endcap_ =
0729         std::make_unique<CLHEP::RandGaussQ>(eng, theThresholdInE_Endcap_, theThresholdSmearing_Endcap_);
0730     smearedThreshold_Barrel_ =
0731         std::make_unique<CLHEP::RandGaussQ>(eng, theThresholdInE_Barrel_, theThresholdSmearing_Barrel_);
0732   }
0733   rengine_ = &eng;
0734   _signal.clear();
0735 }
0736 
0737 // =======================================================================================
0738 //
0739 // Set the drift direction accoring to the Bfield in local det-unit frame
0740 // Works for both barrel and forward pixels.
0741 // Replace the sign convention to fit M.Swartz's formulaes.
0742 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
0743 // parameter value
0744 
0745 LocalVector Phase2TrackerDigitizerAlgorithm::DriftDirection(const Phase2TrackerGeomDetUnit* pixdet,
0746                                                             const GlobalVector& bfield,
0747                                                             const DetId& detId) const {
0748   Frame detFrame(pixdet->surface().position(), pixdet->surface().rotation());
0749   LocalVector Bfield = detFrame.toLocal(bfield);
0750 
0751   float dir_x = 0.0;
0752   float dir_y = 0.0;
0753   float dir_z = 0.0;
0754   float scale = 1.0;
0755 
0756   uint32_t detID = pixdet->geographicalId().rawId();
0757   uint32_t Sub_detid = DetId(detID).subdetId();
0758 
0759   // Read Lorentz angle from DB:
0760   if (use_LorentzAngle_DB_) {
0761     bool isPixel = (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == PixelSubdetector::PixelEndcap);
0762 
0763     float lorentzAngle =
0764         isPixel ? siPixelLorentzAngle_->getLorentzAngle(detId) : siPhase2OTLorentzAngle_->getLorentzAngle(detId);
0765     float alpha2 = std::pow(lorentzAngle, 2);
0766 
0767     dir_x = -(lorentzAngle * Bfield.y() + alpha2 * Bfield.z() * Bfield.x());
0768     dir_y = +(lorentzAngle * Bfield.x() - alpha2 * Bfield.z() * Bfield.y());
0769     dir_z = -(1 + alpha2 * std::pow(Bfield.z(), 2));
0770     scale = (1 + alpha2 * std::pow(Bfield.z(), 2));
0771   } else {
0772     // Read Lorentz angle from cfg file:
0773     float alpha2_Endcap = 0.0;
0774     float alpha2_Barrel = 0.0;
0775     if (alpha2Order_) {
0776       alpha2_Endcap = std::pow(tanLorentzAnglePerTesla_Endcap_, 2);
0777       alpha2_Barrel = std::pow(tanLorentzAnglePerTesla_Barrel_, 2);
0778     }
0779 
0780     if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) {  // barrel layers
0781       dir_x = -(tanLorentzAnglePerTesla_Barrel_ * Bfield.y() + alpha2_Barrel * Bfield.z() * Bfield.x());
0782       dir_y = +(tanLorentzAnglePerTesla_Barrel_ * Bfield.x() - alpha2_Barrel * Bfield.z() * Bfield.y());
0783       dir_z = -(1 + alpha2_Barrel * std::pow(Bfield.z(), 2));
0784       scale = (1 + alpha2_Barrel * std::pow(Bfield.z(), 2));
0785 
0786     } else {  // forward disks
0787       dir_x = -(tanLorentzAnglePerTesla_Endcap_ * Bfield.y() + alpha2_Endcap * Bfield.z() * Bfield.x());
0788       dir_y = +(tanLorentzAnglePerTesla_Endcap_ * Bfield.x() - alpha2_Endcap * Bfield.z() * Bfield.y());
0789       dir_z = -(1 + alpha2_Endcap * std::pow(Bfield.z(), 2));
0790       scale = (1 + alpha2_Endcap * std::pow(Bfield.z(), 2));
0791     }
0792   }
0793 
0794   LocalVector theDriftDirection = LocalVector(dir_x / scale, dir_y / scale, dir_z / scale);
0795   LogDebug("Phase2TrackerDigitizerAlgorithm") << " The drift direction in local coordinate is " << theDriftDirection;
0796   return theDriftDirection;
0797 }
0798 
0799 // =============================================================================
0800 
0801 void Phase2TrackerDigitizerAlgorithm::pixel_inefficiency_db(uint32_t detID) {
0802   signal_map_type& theSignal = _signal[detID];  // check validity
0803 
0804   // Loop over hit pixels, amplitude in electrons, channel = coded row,col
0805   for (auto& s : theSignal) {
0806     std::pair<int, int> ip;
0807     if (pixelFlag_)
0808       ip = PixelDigi::channelToPixel(s.first);  //get pixel pos
0809     else
0810       ip = Phase2TrackerDigi::channelToPixel(s.first);  //get pixel pos
0811 
0812     int row = ip.first;   // X in row
0813     int col = ip.second;  // Y is in col
0814     // transform to ROC index coordinates
0815     if (theSiPixelGainCalibrationService_->isDead(detID, col, row))
0816       s.second.set(0.);  // reset amplitude
0817   }
0818 }
0819 
0820 // ==========================================================================
0821 
0822 void Phase2TrackerDigitizerAlgorithm::module_killing_conf(uint32_t detID) {
0823   bool isbad = false;
0824   int detid = detID;
0825   std::string Module;
0826   for (auto const& det_m : deadModules_) {
0827     int Dead_detID = det_m.getParameter<int>("Dead_detID");
0828     Module = det_m.getParameter<std::string>("Module");
0829     if (detid == Dead_detID) {
0830       isbad = true;
0831       break;
0832     }
0833   }
0834 
0835   if (!isbad)
0836     return;
0837 
0838   signal_map_type& theSignal = _signal[detID];  // check validity
0839   for (auto& s : theSignal) {
0840     std::pair<int, int> ip;
0841     if (pixelFlag_)
0842       ip = PixelDigi::channelToPixel(s.first);
0843     else
0844       ip = Phase2TrackerDigi::channelToPixel(s.first);  //get pixel pos
0845 
0846     if (Module == "whole" || (Module == "tbmA" && ip.first >= 80 && ip.first <= 159) ||
0847         (Module == "tbmB" && ip.first <= 79))
0848       s.second.set(0.);
0849   }
0850 }
0851 // For premixing
0852 void Phase2TrackerDigitizerAlgorithm::loadAccumulator(uint32_t detId, const std::map<int, float>& accumulator) {
0853   auto& theSignal = _signal[detId];
0854   // the input channel is always with PixelDigi definition
0855   // if needed, that has to be converted to Phase2TrackerDigi convention
0856   for (const auto& elem : accumulator) {
0857     auto inserted = theSignal.emplace(elem.first, DigitizerUtility::Amplitude(elem.second, nullptr));
0858     if (!inserted.second) {
0859       throw cms::Exception("LogicError") << "Signal was already set for DetId " << detId;
0860     }
0861   }
0862 }
0863 
0864 void Phase2TrackerDigitizerAlgorithm::digitize(const Phase2TrackerGeomDetUnit* pixdet,
0865                                                std::map<int, DigitizerUtility::DigiSimInfo>& digi_map,
0866                                                const TrackerTopology* tTopo) {
0867   uint32_t detID = pixdet->geographicalId().rawId();
0868   auto it = _signal.find(detID);
0869   if (it == _signal.end())
0870     return;
0871 
0872   const signal_map_type& theSignal = _signal[detID];
0873 
0874   uint32_t Sub_detid = DetId(detID).subdetId();
0875 
0876   float theThresholdInE = 0.;
0877   float theHIPThresholdInE = 0.;
0878   // Define Threshold
0879   if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) {  // Barrel modules
0880     theThresholdInE = addThresholdSmearing_ ? smearedThreshold_Barrel_->fire()             // gaussian smearing
0881                                             : theThresholdInE_Barrel_;                     // no smearing
0882     theHIPThresholdInE = theHIPThresholdInE_Barrel_;
0883   } else {                                                                      // Forward disks modules
0884     theThresholdInE = addThresholdSmearing_ ? smearedThreshold_Endcap_->fire()  // gaussian smearing
0885                                             : theThresholdInE_Endcap_;          // no smearing
0886     theHIPThresholdInE = theHIPThresholdInE_Endcap_;
0887   }
0888 
0889   //  if (addNoise) add_noise(pixdet, theThresholdInE/theNoiseInElectrons_);  // generate noise
0890   if (addNoise_)
0891     add_noise(pixdet);  // generate noise
0892   if (addXtalk_)
0893     add_cross_talk(pixdet);
0894   if (addNoisyPixels_)
0895     add_noisy_cells(pixdet, theHIPThresholdInE / theElectronPerADC_);
0896 
0897   // Do only if needed
0898   if (addPixelInefficiency_ && !theSignal.empty()) {
0899     if (use_ineff_from_db_)
0900       pixel_inefficiency_db(detID);
0901     else
0902       pixel_inefficiency(subdetEfficiencies_, pixdet, tTopo);
0903   }
0904   if (use_module_killing_) {
0905     if (use_deadmodule_DB_)  // remove dead modules using DB
0906       module_killing_DB(pixdet);
0907     else  // remove dead modules using the list in cfg file
0908       module_killing_conf(detID);
0909   }
0910 
0911   // Digitize if the signal is greater than threshold
0912   for (auto const& s : theSignal) {
0913     const DigitizerUtility::Amplitude& sig_data = s.second;
0914     float signalInElectrons = sig_data.ampl();
0915 
0916     const auto& info_list = sig_data.simInfoList();
0917     const DigitizerUtility::SimHitInfo* hitInfo = nullptr;
0918     if (!info_list.empty())
0919       hitInfo = std::max_element(info_list.begin(), info_list.end())->second.get();
0920 
0921     if (isAboveThreshold(hitInfo, signalInElectrons, theThresholdInE)) {  // check threshold
0922       DigitizerUtility::DigiSimInfo info;
0923       info.sig_tot = convertSignalToAdc(detID, signalInElectrons, theThresholdInE);  // adc
0924       info.ot_bit = signalInElectrons > theHIPThresholdInE ? true : false;
0925       if (makeDigiSimLinks_) {
0926         for (auto const& l : sig_data.simInfoList()) {
0927           float charge_frac = l.first / signalInElectrons;
0928           if (l.first > -5.0)
0929             info.simInfoList.push_back({charge_frac, l.second.get()});
0930         }
0931       }
0932       digi_map.insert({s.first, info});
0933     }
0934   }
0935 }
0936 //
0937 // Scale the Signal using Dual Slope option
0938 //
0939 int Phase2TrackerDigitizerAlgorithm::convertSignalToAdc(uint32_t detID, float signal_in_elec, float threshold) {
0940   int signal_in_adc;
0941   int temp_signal;
0942   const int max_limit = 10;
0943   if (thePhase2ReadoutMode_ == 0) {
0944     signal_in_adc = theAdcFullScale_;
0945   } else {
0946     if (thePhase2ReadoutMode_ == -1) {
0947       temp_signal = std::min(static_cast<int>(signal_in_elec / theElectronPerADC_), theAdcFullScale_);
0948     } else {
0949       // calculate the kink point and the slope
0950       int dualslope_param = std::min(std::abs(thePhase2ReadoutMode_), max_limit);
0951       int kink_point = static_cast<int>(theAdcFullScale_ / 2) + 1;
0952       // C-ROC: first valid ToT code above threshold is 0b0000
0953       temp_signal = std::floor((signal_in_elec - threshold) / theElectronPerADC_);
0954       if (temp_signal > kink_point)
0955         temp_signal = std::floor((temp_signal - kink_point) / (pow(2, dualslope_param - 1))) + kink_point;
0956     }
0957     signal_in_adc = std::min(temp_signal, theAdcFullScale_);
0958     LogTrace("Phase2TrackerDigitizerAlgorithm")
0959         << " DetId " << detID << " signal_in_elec " << signal_in_elec << " threshold " << threshold
0960         << " signal_above_thr " << signal_in_elec - threshold << " temp conversion "
0961         << std::floor((signal_in_elec - threshold) / theElectronPerADC_) + 1 << " signal after slope correction "
0962         << temp_signal << " signal_in_adc " << signal_in_adc;
0963   }
0964   return signal_in_adc;
0965 }
0966 float Phase2TrackerDigitizerAlgorithm::calcQ(float x) {
0967   constexpr float p1 = 12.5f;
0968   constexpr float p2 = 0.2733f;
0969   constexpr float p3 = 0.147f;
0970 
0971   auto xx = std::min(0.5f * x * x, p1);
0972   return 0.5f * (1.f - std::copysign(std::sqrt(1.f - unsafe_expf<4>(-xx * (1.f + p2 / (1.f + p3 * xx)))), x));
0973 }