CMS 3D CMS Logo

Phase2TrackerDigitizerAlgorithm.cc
Go to the documentation of this file.
1 #include <cmath>
2 #include <iostream>
3 #include <typeinfo>
4 
5 #include "CLHEP/Random/RandGaussQ.h"
26 
27 using namespace edm;
28 
29 namespace ph2tkdigialgo {
30  // Mass in MeV
31  constexpr double m_pion = 139.571;
32  constexpr double m_kaon = 493.677;
33  constexpr double m_electron = 0.511;
34  constexpr double m_muon = 105.658;
35  constexpr double m_proton = 938.272;
36 } // namespace ph2tkdigialgo
37 
39  const edm::ParameterSet& conf_specific,
41  : _signal(),
42  makeDigiSimLinks_(conf_common.getUntrackedParameter<bool>("makeDigiSimLinks", true)),
43  use_ineff_from_db_(conf_specific.getParameter<bool>("Inefficiency_DB")),
44  use_module_killing_(conf_specific.getParameter<bool>("KillModules")), // boolean to kill or not modules
45  use_deadmodule_DB_(conf_specific.getParameter<bool>("DeadModules_DB")), // boolean to access dead modules from DB
46  // boolean to access Lorentz angle from DB
47  use_LorentzAngle_DB_(conf_specific.getParameter<bool>("LorentzAngle_DB")),
48 
49  // get dead module from cfg file
50  deadModules_(use_deadmodule_DB_ ? Parameters() : conf_specific.getParameter<Parameters>("DeadModules")),
51 
52  // Common pixel parameters
53  // These are parameters which are not likely to be changed
54  GeVperElectron_(3.61E-09), // 1 electron(3.61eV, 1keV(277e, mod 9/06 d.k.
55  alpha2Order_(conf_specific.getParameter<bool>("Alpha2Order")), // switch on/off of E.B effect
56  addXtalk_(conf_specific.getParameter<bool>("AddXTalk")),
57  // Interstrip Coupling - Not used in PixelDigitizerAlgorithm
58  interstripCoupling_(conf_specific.getParameter<double>("InterstripCoupling")),
59 
60  Sigma0_(conf_specific.getParameter<double>("SigmaZero")), // Charge diffusion constant 7->3.7
61  SigmaCoeff_(conf_specific.getParameter<double>("SigmaCoeff")), // delta in the diffusion across the strip pitch
62  // (set between 0 to 0.9, 0-->flat Sigma0, 1-->Sigma_min=0 & Sigma_max=2*Sigma0
63  // D.B.: Dist300 replaced by moduleThickness, may not work with partially depleted sensors but works otherwise
64  // Dist300(0.0300), // normalized to 300micron Silicon
65 
66  // Charge integration spread on the collection plane
67  clusterWidth_(conf_specific.getParameter<double>("ClusterWidth")),
68 
69  // Allowed modes of readout which has following values :
70  // 0 ---> Digital or binary readout
71  // -1 ---> Analog readout, current digitizer (Inner Pixel) (TDR version) with no threshold subtraction
72  // Analog readout with dual slope with the "second" slope being 1/2^(n-1) and threshold subtraction (n = 1, 2, 3,4)
73  thePhase2ReadoutMode_(conf_specific.getParameter<int>("Phase2ReadoutMode")),
74 
75  // ADC calibration 1adc count(135e.
76  // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
77  // Be careful, this parameter is also used in SiPixelDet.cc to
78  // calculate the noise in adc counts from noise in electrons.
79  // Both defaults should be the same.
80  theElectronPerADC_(conf_specific.getParameter<double>("ElectronPerAdc")),
81 
82  // ADC saturation value, 255(8bit adc.
83  theAdcFullScale_(conf_specific.getParameter<int>("AdcFullScale")),
84 
85  // Noise in electrons:
86  // Pixel cell noise, relevant for generating noisy pixels
87  theNoiseInElectrons_(conf_specific.getParameter<double>("NoiseInElectrons")),
88 
89  // Fill readout noise, including all readout chain, relevant for smearing
90  theReadoutNoise_(conf_specific.getParameter<double>("ReadoutNoiseInElec")),
91 
92  // Threshold in units of noise:
93  // thePixelThreshold(conf.getParameter<double>("ThresholdInNoiseUnits")),
94  // Pixel threshold in electron units.
95  theThresholdInE_Endcap_(conf_specific.getParameter<double>("ThresholdInElectrons_Endcap")),
96  theThresholdInE_Barrel_(conf_specific.getParameter<double>("ThresholdInElectrons_Barrel")),
97 
98  // Add threshold gaussian smearing:
99  theThresholdSmearing_Endcap_(conf_specific.getParameter<double>("ThresholdSmearing_Endcap")),
100  theThresholdSmearing_Barrel_(conf_specific.getParameter<double>("ThresholdSmearing_Barrel")),
101 
102  // Add HIP Threshold in electron units.
103  theHIPThresholdInE_Endcap_(conf_specific.getParameter<double>("HIPThresholdInElectrons_Endcap")),
104  theHIPThresholdInE_Barrel_(conf_specific.getParameter<double>("HIPThresholdInElectrons_Barrel")),
105 
106  // theTofCut 12.5, cut in particle TOD +/- 12.5ns
107  theTofLowerCut_(conf_specific.getParameter<double>("TofLowerCut")),
108  theTofUpperCut_(conf_specific.getParameter<double>("TofUpperCut")),
109 
110  // Get the Lorentz angle from the cfg file:
111  tanLorentzAnglePerTesla_Endcap_(
112  use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Endcap")),
113  tanLorentzAnglePerTesla_Barrel_(
114  use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Barrel")),
115 
116  // Add noise
117  addNoise_(conf_specific.getParameter<bool>("AddNoise")),
118 
119  // Add noisy pixels
120  addNoisyPixels_(conf_specific.getParameter<bool>("AddNoisyPixels")),
121 
122  // Fluctuate charge in track subsegments
123  fluctuateCharge_(conf_specific.getUntrackedParameter<bool>("FluctuateCharge", true)),
124 
125  // Control the pixel inefficiency
126  addPixelInefficiency_(conf_specific.getParameter<bool>("AddInefficiency")),
127 
128  // Add threshold gaussian smearing:
129  addThresholdSmearing_(conf_specific.getParameter<bool>("AddThresholdSmearing")),
130 
131  // Add some pseudo-red damage
132  pseudoRadDamage_(conf_specific.exists("PseudoRadDamage") ? conf_specific.getParameter<double>("PseudoRadDamage")
133  : double(0.0)),
134  pseudoRadDamageRadius_(conf_specific.exists("PseudoRadDamageRadius")
135  ? conf_specific.getParameter<double>("PseudoRadDamageRadius")
136  : double(0.0)),
137  // Add pixel radiation damage
138  useChargeReweighting_(conf_specific.getParameter<bool>("UseReweighting")),
139  theSiPixelChargeReweightingAlgorithm_(
140  useChargeReweighting_ ? std::make_unique<SiPixelChargeReweightingAlgorithm>(conf_specific, iC) : nullptr),
141 
142  // delta cutoff in MeV, has to be same as in OSCAR(0.030/cmsim=1.0 MeV
143  // tMax(0.030), // In MeV.
144  // tMax(conf.getUntrackedParameter<double>("DeltaProductionCut",0.030)),
145  tMax_(conf_common.getParameter<double>("DeltaProductionCut")),
146 
147  badPixels_(conf_specific.getParameter<Parameters>("CellsToKill")),
148 
149  fluctuate_(fluctuateCharge_ ? std::make_unique<SiG4UniversalFluctuation>() : nullptr),
150  theNoiser_(addNoise_ ? std::make_unique<GaussianTailNoiseGenerator>() : nullptr),
151  theSiPixelGainCalibrationService_(
152  use_ineff_from_db_ ? std::make_unique<SiPixelGainCalibrationOfflineSimService>(conf_specific, iC) : nullptr),
153  subdetEfficiencies_(conf_specific) {
154  LogDebug("Phase2TrackerDigitizerAlgorithm")
155  << "Phase2TrackerDigitizerAlgorithm constructed\n"
156  << "Configuration parameters:\n"
157  << "Threshold/Gain = "
158  << "threshold in electron Endcap = " << theThresholdInE_Endcap_
159  << "\nthreshold in electron Barrel = " << theThresholdInE_Barrel_ << " ElectronPerADC " << theElectronPerADC_
160  << " ADC Scale (in bits) " << theAdcFullScale_ << " The delta cut-off is set to " << tMax_ << " pix-inefficiency "
162 }
163 
165  LogDebug("Phase2TrackerDigitizerAlgorithm") << "Phase2TrackerDigitizerAlgorithm deleted";
166 }
167 
169  barrel_efficiencies = conf.getParameter<std::vector<double> >("EfficiencyFactors_Barrel");
170  endcap_efficiencies = conf.getParameter<std::vector<double> >("EfficiencyFactors_Endcap");
171 }
172 
173 void Phase2TrackerDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
174  std::vector<PSimHit>::const_iterator inputEnd,
175  const size_t inputBeginGlobalIndex,
176  const uint32_t tofBin,
177  const Phase2TrackerGeomDetUnit* pixdet,
178  const GlobalVector& bfield) {
179  // produce SignalPoint's for all SimHit's in detector
180  // Loop over hits
181  uint32_t detId = pixdet->geographicalId().rawId();
182  size_t simHitGlobalIndex = inputBeginGlobalIndex; // This needs to be stored to create the digi-sim link later
183 
184  // find the relevant hits
185  std::vector<PSimHit> matchedSimHits;
186  std::copy_if(inputBegin, inputEnd, std::back_inserter(matchedSimHits), [detId](auto const& hit) -> bool {
187  return hit.detUnitId() == detId;
188  });
189  // loop over a much reduced set of SimHits
190  for (auto const& hit : matchedSimHits) {
191  LogDebug("Phase2DigitizerAlgorithm") << hit.particleType() << " " << hit.pabs() << " " << hit.energyLoss() << " "
192  << hit.tof() << " " << hit.trackId() << " " << hit.processType() << " "
193  << hit.detUnitId() << hit.entryPoint() << " " << hit.exitPoint();
194  double signalScale = 1.0;
195  // fill collection_points for this SimHit, indpendent of topology
196  if (select_hit(hit, (pixdet->surface().toGlobal(hit.localPosition()).mag() * c_inv), signalScale)) {
197  const auto& ionization_points = primary_ionization(hit); // fills ionization_points
198 
199  // transforms ionization_points -> collection_points
200  const auto& collection_points = drift(hit, pixdet, bfield, ionization_points);
201 
202  // compute induced signal on readout elements and add to _signal
203  // hit needed only for SimHit<-->Digi link
204  induce_signal(inputBegin, hit, simHitGlobalIndex, inputBeginGlobalIndex, tofBin, pixdet, collection_points);
205  LogDebug("Phase2DigitizerAlgorithm") << "Induced signal was computed";
206  }
207  ++simHitGlobalIndex;
208  }
209 }
210 // =================================================================
211 //
212 // Generate primary ionization along the track segment.
213 // Divide the track into small sub-segments
214 //
215 // =================================================================
216 std::vector<digitizerUtility::EnergyDepositUnit> Phase2TrackerDigitizerAlgorithm::primary_ionization(
217  const PSimHit& hit) const {
218  // Straight line approximation for trajectory inside active media
219  constexpr float SegmentLength = 0.0010; // in cm (10 microns)
220  // Get the 3D segment direction vector
221  LocalVector direction = hit.exitPoint() - hit.entryPoint();
222 
223  float eLoss = hit.energyLoss(); // Eloss in GeV
224  float length = direction.mag(); // Track length in Silicon
225 
226  int NumberOfSegments = static_cast<int>(length / SegmentLength); // Number of segments
227  if (NumberOfSegments < 1)
228  NumberOfSegments = 1;
229  LogDebug("Phase2TrackerDigitizerAlgorithm")
230  << "enter primary_ionzation " << NumberOfSegments << " shift = " << hit.exitPoint().x() - hit.entryPoint().x()
231  << " " << hit.exitPoint().y() - hit.entryPoint().y() << " " << hit.exitPoint().z() - hit.entryPoint().z() << " "
232  << hit.particleType() << " " << hit.pabs();
233 
234  std::vector<float> elossVector;
235  elossVector.reserve(NumberOfSegments);
236  if (fluctuateCharge_) {
237  // Generate fluctuated charge points
238  elossVector = fluctuateEloss(hit.particleType(), hit.pabs(), eLoss, length, NumberOfSegments);
239  } else {
240  float averageEloss = eLoss / NumberOfSegments;
241  elossVector.resize(NumberOfSegments, averageEloss);
242  }
243 
244  std::vector<digitizerUtility::EnergyDepositUnit> ionization_points;
245  ionization_points.reserve(NumberOfSegments); // set size
246  // loop over segments
247  for (size_t i = 0; i < elossVector.size(); ++i) {
248  // Divide the segment into equal length subsegments
249  Local3DPoint point = hit.entryPoint() + ((i + 0.5) / NumberOfSegments) * direction;
250  float energy = elossVector[i] / GeVperElectron_; // Convert charge to elec.
251 
252  digitizerUtility::EnergyDepositUnit edu(energy, point); // define position,energy point
253  ionization_points.push_back(edu); // save
254  LogDebug("Phase2TrackerDigitizerAlgorithm")
255  << "For index = " << i << " EnergyDepositUnit-x = " << edu.x() << " EnergyDepositUnit-y = " << edu.y()
256  << " EnergyDepositUnit-z = " << edu.z() << " EnergyDepositUnit-energy = " << edu.energy();
257  }
258  return ionization_points;
259 }
260 //==============================================================================
261 //
262 // Fluctuate the charge comming from a small (10um) track segment.
263 // Use the G4 routine. For mip pions for the moment.
264 //
265 //==============================================================================
267  int pid, float particleMomentum, float eloss, float length, int NumberOfSegs) const {
268  double particleMass = ph2tkdigialgo::m_pion; // Mass in MeV, assume pion
269  pid = std::abs(pid);
270  if (pid != 211) { // Mass in MeV
271  if (pid == 11)
272  particleMass = ph2tkdigialgo::m_electron;
273  else if (pid == 13)
274  particleMass = ph2tkdigialgo::m_muon;
275  else if (pid == 321)
276  particleMass = ph2tkdigialgo::m_kaon;
277  else if (pid == 2212)
278  particleMass = ph2tkdigialgo::m_proton;
279  }
280  // What is the track segment length.
281  float segmentLength = length / NumberOfSegs;
282 
283  // Generate charge fluctuations.
284  float sum = 0.;
285  double segmentEloss = (1000. * eloss) / NumberOfSegs; //eloss in MeV
286  std::vector<float> elossVector;
287  elossVector.reserve(NumberOfSegs);
288  for (int i = 0; i < NumberOfSegs; ++i) {
289  // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
290  // track segment length in mm, segment eloss in MeV
291  // Returns fluctuated eloss in MeV
292  double deltaCutoff = tMax_; // the cutoff is sometimes redefined inside, so fix it.
293  float de = fluctuate_->SampleFluctuations(particleMomentum * 1000.,
294  particleMass,
295  deltaCutoff,
296  segmentLength * 10.,
297  segmentEloss,
298  rengine_) /
299  1000.; //convert to GeV
300  elossVector.push_back(de);
301  sum += de;
302  }
303  if (sum > 0.) { // if fluctuations give eloss>0.
304  // Rescale to the same total eloss
305  float ratio = eloss / sum;
307  std::begin(elossVector), std::end(elossVector), std::begin(elossVector), [&ratio](auto const& c) -> float {
308  return c * ratio;
309  }); // use a simple lambda expression
310  } else { // if fluctuations gives 0 eloss
311  float averageEloss = eloss / NumberOfSegs;
312  elossVector.resize(NumberOfSegs, averageEloss);
313  }
314  return elossVector;
315 }
316 
317 // ======================================================================
318 //
319 // Drift the charge segments to the sensor surface (collection plane)
320 // Include the effect of E-field and B-field
321 //
322 // =====================================================================
323 std::vector<digitizerUtility::SignalPoint> Phase2TrackerDigitizerAlgorithm::drift(
324  const PSimHit& hit,
325  const Phase2TrackerGeomDetUnit* pixdet,
326  const GlobalVector& bfield,
327  const std::vector<digitizerUtility::EnergyDepositUnit>& ionization_points) const {
328  LogDebug("Phase2TrackerDigitizerAlgorithm") << "Enter drift ";
329 
330  std::vector<digitizerUtility::SignalPoint> collection_points;
331  collection_points.reserve(ionization_points.size()); // set size
332  LocalVector driftDir = driftDirection(pixdet, bfield, hit.detUnitId()); // get the charge drift direction
333  if (driftDir.z() == 0.) {
334  LogWarning("Phase2TrackerDigitizerAlgorithm") << " pxlx: drift in z is zero ";
335  return collection_points;
336  }
337 
338  float TanLorenzAngleX = driftDir.x(); // tangent of Lorentz angle
339  float TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
340  float dir_z = driftDir.z(); // The z drift direction
341  float CosLorenzAngleX = 1. / std::sqrt(1. + std::pow(TanLorenzAngleX, 2)); // cosine to estimate the path length
342  float CosLorenzAngleY = 1.;
343  if (alpha2Order_) {
344  TanLorenzAngleY = driftDir.y();
345  CosLorenzAngleY = 1. / std::sqrt(1. + std::pow(TanLorenzAngleY, 2)); // cosine
346  }
347 
348  float moduleThickness = pixdet->specificSurface().bounds().thickness();
349  float stripPitch = pixdet->specificTopology().pitch().first;
350 
351  LogDebug("Phase2TrackerDigitizerAlgorithm")
352  << " Lorentz Tan-X " << TanLorenzAngleX << " Lorentz Tan-Y " << TanLorenzAngleY << " Lorentz Cos-X "
353  << CosLorenzAngleX << " Lorentz Cos-Y " << CosLorenzAngleY
354  << " ticknes * Lorentz Tan-X = " << moduleThickness * TanLorenzAngleX << " drift direction " << driftDir;
355 
356  for (auto const& val : ionization_points) {
357  // position
358  float SegX = val.x(), SegY = val.y(), SegZ = val.z();
359 
360  // Distance from the collection plane
361  // DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
362  // Include explixitely the E drift direction (for CMS dir_z=-1)
363 
364  // Distance between charge generation and collection
365  float driftDistance = moduleThickness / 2. - (dir_z * SegZ); // Drift to -z
366 
367  if (driftDistance < 0.)
368  driftDistance = 0.;
369  else if (driftDistance > moduleThickness)
370  driftDistance = moduleThickness;
371 
372  // Assume full depletion now, partial depletion will come later.
373  float XDriftDueToMagField = driftDistance * TanLorenzAngleX;
374  float YDriftDueToMagField = driftDistance * TanLorenzAngleY;
375 
376  // Shift cloud center
377  float CloudCenterX = SegX + XDriftDueToMagField;
378  float CloudCenterY = SegY + YDriftDueToMagField;
379 
380  // Calculate how long is the charge drift path
381  // Actual Drift Lentgh
382  float driftLength =
383  std::sqrt(std::pow(driftDistance, 2) + std::pow(XDriftDueToMagField, 2) + std::pow(YDriftDueToMagField, 2));
384 
385  // What is the charge diffusion after this path
386  // Sigma0=0.00037 is for 300um thickness (make sure moduleThickness is in [cm])
387  float Sigma = std::sqrt(driftLength / moduleThickness) * Sigma0_ * moduleThickness / 0.0300;
388  // D.B.: sigmaCoeff=0 means no modulation
389  if (SigmaCoeff_)
390  Sigma *= (SigmaCoeff_ * std::pow(cos(SegX * M_PI / stripPitch), 2) + 1);
391  // NB: divided by 4 to get a periodicity of stripPitch
392 
393  // Project the diffusion sigma on the collection plane
394  float Sigma_x = Sigma / CosLorenzAngleX;
395  float Sigma_y = Sigma / CosLorenzAngleY;
396 
397  // Insert a charge loss due to Rad Damage here
398  float energyOnCollector = val.energy(); // The energy that reaches the collector
399 
400  // pseudoRadDamage
401  if (pseudoRadDamage_) {
402  float moduleRadius = pixdet->surface().position().perp();
403  if (moduleRadius <= pseudoRadDamageRadius_) {
404  float kValue = pseudoRadDamage_ / std::pow(moduleRadius, 2);
405  energyOnCollector *= exp(-1 * kValue * driftDistance / moduleThickness);
406  }
407  }
408  LogDebug("Phase2TrackerDigitizerAlgorithm")
409  << "Dift DistanceZ = " << driftDistance << " module thickness = " << moduleThickness
410  << " Start Energy = " << val.energy() << " Energy after loss= " << energyOnCollector;
411  digitizerUtility::SignalPoint sp(CloudCenterX, CloudCenterY, Sigma_x, Sigma_y, hit.tof(), energyOnCollector);
412 
413  // Load the Charge distribution parameters
414  collection_points.push_back(sp);
415  }
416  return collection_points;
417 }
418 
419 // ====================================================================
420 //
421 // Induce the signal on the collection plane of the active sensor area.
423  std::vector<PSimHit>::const_iterator inputBegin,
424  const PSimHit& hit,
425  const size_t hitIndex,
426  const size_t firstHitIndex,
427  const uint32_t tofBin,
428  const Phase2TrackerGeomDetUnit* pixdet,
429  const std::vector<digitizerUtility::SignalPoint>& collection_points) {
430  // X - Rows, Left-Right, 160, (1.6cm) for barrel
431  // Y - Columns, Down-Up, 416, (6.4cm)
432  const Phase2TrackerTopology* topol = &pixdet->specificTopology();
433  uint32_t detID = pixdet->geographicalId().rawId();
434  signal_map_type& theSignal = _signal[detID];
435 
436  LogDebug("Phase2TrackerDigitizerAlgorithm")
437  << "Enter induce_signal, Pitch size is " << topol->pitch().first << " cm vs " << topol->pitch().second << " cm";
438 
439  // local map to store pixels hit by 1 Hit.
440  using hit_map_type = std::map<int, float, std::less<int> >;
441  hit_map_type hit_signal;
442 
443  // Assign signals to readout channels and store sorted by channel number
444  // Iterate over collection points on the collection plane
445  for (auto const& v : collection_points) {
446  float CloudCenterX = v.position().x(); // Charge position in x
447  float CloudCenterY = v.position().y(); // in y
448  float SigmaX = v.sigma_x(); // Charge spread in x
449  float SigmaY = v.sigma_y(); // in y
450  float Charge = v.amplitude(); // Charge amplitude
451 
452  LogDebug("Phase2TrackerDigitizerAlgorithm") << " cloud " << v.position().x() << " " << v.position().y() << " "
453  << v.sigma_x() << " " << v.sigma_y() << " " << v.amplitude();
454 
455  // Find the maximum cloud spread in 2D plane , assume 3*sigma
456  float CloudRight = CloudCenterX + clusterWidth_ * SigmaX;
457  float CloudLeft = CloudCenterX - clusterWidth_ * SigmaX;
458  float CloudUp = CloudCenterY + clusterWidth_ * SigmaY;
459  float CloudDown = CloudCenterY - clusterWidth_ * SigmaY;
460 
461  // Define 2D cloud limit points
462  LocalPoint PointRightUp = LocalPoint(CloudRight, CloudUp);
463  LocalPoint PointLeftDown = LocalPoint(CloudLeft, CloudDown);
464 
465  // This points can be located outside the sensor area.
466  // The conversion to measurement point does not check for that
467  // so the returned pixel index might be wrong (outside range).
468  // We rely on the limits check below to fix this.
469  // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
470  // PIXEL ARE IS LOST, it should not be collected.
471 
472  // Convert the 2D points to pixel indices
473  MeasurementPoint mp = topol->measurementPosition(PointRightUp);
474  int IPixRightUpX = static_cast<int>(std::floor(mp.x())); // cast reqd.
475  int IPixRightUpY = static_cast<int>(std::floor(mp.y()));
476  LogDebug("Phase2TrackerDigitizerAlgorithm")
477  << " right-up " << PointRightUp << " " << mp.x() << " " << mp.y() << " " << IPixRightUpX << " " << IPixRightUpY;
478 
479  mp = topol->measurementPosition(PointLeftDown);
480  int IPixLeftDownX = static_cast<int>(std::floor(mp.x()));
481  int IPixLeftDownY = static_cast<int>(std::floor(mp.y()));
482  LogDebug("Phase2TrackerDigitizerAlgorithm") << " left-down " << PointLeftDown << " " << mp.x() << " " << mp.y()
483  << " " << IPixLeftDownX << " " << IPixLeftDownY;
484 
485  // Check detector limits to correct for pixels outside range.
486  int numColumns = topol->ncolumns(); // det module number of cols&rows
487  int numRows = topol->nrows();
488 
489  IPixRightUpX = numRows > IPixRightUpX ? IPixRightUpX : numRows - 1;
490  IPixRightUpY = numColumns > IPixRightUpY ? IPixRightUpY : numColumns - 1;
491  IPixLeftDownX = 0 < IPixLeftDownX ? IPixLeftDownX : 0;
492  IPixLeftDownY = 0 < IPixLeftDownY ? IPixLeftDownY : 0;
493 
494  // First integrate charge strips in x
495  hit_map_type x;
496  for (int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) { // loop over x index
497  float xLB, LowerBound;
498  // Why is set to 0 if ix=0, does it meen that we accept charge
499  // outside the sensor?
500  if (ix == 0 || SigmaX == 0.) { // skip for surface segemnts
501  LowerBound = 0.;
502  } else {
503  mp = MeasurementPoint(ix, 0.0);
504  xLB = topol->localPosition(mp).x();
505  LowerBound = 1 - calcQ((xLB - CloudCenterX) / SigmaX);
506  }
507 
508  float xUB, UpperBound;
509  if (ix == numRows - 1 || SigmaX == 0.) {
510  UpperBound = 1.;
511  } else {
512  mp = MeasurementPoint(ix + 1, 0.0);
513  xUB = topol->localPosition(mp).x();
514  UpperBound = 1. - calcQ((xUB - CloudCenterX) / SigmaX);
515  }
516  float TotalIntegrationRange = UpperBound - LowerBound; // get strip
517  x.emplace(ix, TotalIntegrationRange); // save strip integral
518  }
519 
520  // Now integrate strips in y
521  hit_map_type y;
522  for (int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) { // loop over y index
523  float yLB, LowerBound;
524  if (iy == 0 || SigmaY == 0.) {
525  LowerBound = 0.;
526  } else {
527  mp = MeasurementPoint(0.0, iy);
528  yLB = topol->localPosition(mp).y();
529  LowerBound = 1. - calcQ((yLB - CloudCenterY) / SigmaY);
530  }
531 
532  float yUB, UpperBound;
533  if (iy == numColumns - 1 || SigmaY == 0.) {
534  UpperBound = 1.;
535  } else {
536  mp = MeasurementPoint(0.0, iy + 1);
537  yUB = topol->localPosition(mp).y();
538  UpperBound = 1. - calcQ((yUB - CloudCenterY) / SigmaY);
539  }
540 
541  float TotalIntegrationRange = UpperBound - LowerBound;
542  y.emplace(iy, TotalIntegrationRange); // save strip integral
543  }
544 
545  // Get the 2D charge integrals by folding x and y strips
546  for (int ix = IPixLeftDownX; ix <= IPixRightUpX; ++ix) { // loop over x index
547  for (int iy = IPixLeftDownY; iy <= IPixRightUpY; ++iy) { // loop over y index
548  float ChargeFraction = Charge * x[ix] * y[iy];
549  int chanFired = -1;
550  if (ChargeFraction > 0.) {
551  chanFired =
553  // Load the amplitude
554  hit_signal[chanFired] += ChargeFraction;
555  }
556 
557  mp = MeasurementPoint(ix, iy);
558  LocalPoint lp = topol->localPosition(mp);
559  int chan = topol->channel(lp);
560 
561  LogDebug("Phase2TrackerDigitizerAlgorithm")
562  << " pixel " << ix << " " << iy << " - "
563  << " " << chanFired << " " << ChargeFraction << " " << mp.x() << " " << mp.y() << " " << lp.x() << " "
564  << lp.y() << " " // givex edge position
565  << chan; // edge belongs to previous ?
566  }
567  }
568  }
569  // Fill the global map with all hit pixels from this event
570  bool reweighted = false;
571  size_t referenceIndex4CR = 0;
572 
573  if (useChargeReweighting_) {
574  if (hit.processType() == 0) {
575  referenceIndex4CR = hitIndex;
576  reweighted =
578  hit_signal,
579  hitIndex,
580  referenceIndex4CR,
581  tofBin,
582  topol,
583  detID,
584  theSignal,
585  hit.processType(),
587  } else {
588  // If it's not the primary particle, use the first hit in the collection as SimHit, which should be the corresponding primary.
589  referenceIndex4CR = firstHitIndex;
590  reweighted =
592  hit_signal,
593  hitIndex,
594  referenceIndex4CR,
595  tofBin,
596  topol,
597  detID,
598  theSignal,
599  hit.processType(),
601  }
602  }
603 
604  float corr_time = hit.tof() - pixdet->surface().toGlobal(hit.localPosition()).mag() * c_inv;
605  if (!reweighted) {
606  for (auto const& hit_s : hit_signal) {
607  int chan = hit_s.first;
609  hit_s.second, &hit, hit_s.second, corr_time, hitIndex, tofBin)
610  : digitizerUtility::Ph2Amplitude(hit_s.second, nullptr, hit_s.second));
611  }
612  }
613 } // end of induce_signal function
614 
615 // ======================================================================
616 //
617 // Add electronic noise to pixel charge
618 //
619 // ======================================================================
621  uint32_t detID = pixdet->geographicalId().rawId();
622  signal_map_type& theSignal = _signal[detID];
623  for (auto& s : theSignal) {
624  float noise = gaussDistribution_->fire();
625  if ((s.second.ampl() + noise) < 0.)
626  s.second.set(0);
627  else
628  s.second += noise;
629  }
630 }
631 
632 // ======================================================================
633 //
634 // Add Cross-talk contribution
635 //
636 // ======================================================================
638  uint32_t detID = pixdet->geographicalId().rawId();
639  signal_map_type& theSignal = _signal[detID];
640  signal_map_type signalNew;
641  const Phase2TrackerTopology* topol = &pixdet->specificTopology();
642  int numRows = topol->nrows();
643 
644  for (auto& s : theSignal) {
645  float signalInElectrons = s.second.ampl(); // signal in electrons
646 
647  std::pair<int, int> hitChan;
648  if (pixelFlag_)
649  hitChan = PixelDigi::channelToPixel(s.first);
650  else
651  hitChan = Phase2TrackerDigi::channelToPixel(s.first);
652 
653  float signalInElectrons_Xtalk = signalInElectrons * interstripCoupling_;
654  // subtract the charge which will be shared
655  s.second.set(signalInElectrons - signalInElectrons_Xtalk);
656 
657  if (hitChan.first != 0) {
658  auto XtalkPrev = std::make_pair(hitChan.first - 1, hitChan.second);
659  int chanXtalkPrev = pixelFlag_ ? PixelDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second)
660  : Phase2TrackerDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second);
661  signalNew.emplace(chanXtalkPrev, digitizerUtility::Ph2Amplitude(signalInElectrons_Xtalk, nullptr, -1.0));
662  }
663  if (hitChan.first < numRows - 1) {
664  auto XtalkNext = std::make_pair(hitChan.first + 1, hitChan.second);
665  int chanXtalkNext = pixelFlag_ ? PixelDigi::pixelToChannel(XtalkNext.first, XtalkNext.second)
666  : Phase2TrackerDigi::pixelToChannel(XtalkNext.first, XtalkNext.second);
667  signalNew.emplace(chanXtalkNext, digitizerUtility::Ph2Amplitude(signalInElectrons_Xtalk, nullptr, -1.0));
668  }
669  }
670  for (auto const& l : signalNew) {
671  int chan = l.first;
672  auto iter = theSignal.find(chan);
673  if (iter != theSignal.end()) {
674  theSignal[chan] += l.second.ampl();
675  } else {
676  theSignal.emplace(chan, digitizerUtility::Ph2Amplitude(l.second.ampl(), nullptr, -1.0));
677  }
678  }
679 }
680 
681 // ======================================================================
682 //
683 // Add noise on non-hit cells
684 //
685 // ======================================================================
687  uint32_t detID = pixdet->geographicalId().rawId();
688  signal_map_type& theSignal = _signal[detID];
689  const Phase2TrackerTopology* topol = &pixdet->specificTopology();
690 
691  int numColumns = topol->ncolumns(); // det module number of cols&rows
692  int numRows = topol->nrows();
693 
694  int numberOfPixels = numRows * numColumns;
695  std::map<int, float, std::less<int> > otherPixels;
696 
697  theNoiser_->generate(numberOfPixels,
698  thePixelThreshold, //thr. in un. of nois
699  theNoiseInElectrons_, // noise in elec.
700  otherPixels,
701  rengine_);
702 
703  LogDebug("Phase2TrackerDigitizerAlgorithm")
704  << " Add noisy pixels " << numRows << " " << numColumns << " " << theNoiseInElectrons_ << " "
705  << theThresholdInE_Endcap_ << " " << theThresholdInE_Barrel_ << " " << numberOfPixels << " "
706  << otherPixels.size();
707 
708  // Add noisy pixels
709  for (auto const& el : otherPixels) {
710  int iy = el.first / numRows;
711  if (iy < 0 || iy > numColumns - 1)
712  LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in iy " << iy;
713 
714  int ix = el.first - iy * numRows;
715  if (ix < 0 || ix > numRows - 1)
716  LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in ix " << ix;
717 
719 
720  LogDebug("Phase2TrackerDigitizerAlgorithm")
721  << " Storing noise = " << el.first << " " << el.second << " " << ix << " " << iy << " " << chan;
722 
723  if (theSignal[chan] == 0)
724  theSignal[chan] = digitizerUtility::Ph2Amplitude(el.second, nullptr, -1.);
725  }
726 }
727 // ============================================================================
728 //
729 // Simulate the readout inefficiencies.
730 // Delete a selected number of single pixels, dcols and rocs.
732  const Phase2TrackerGeomDetUnit* pixdet,
733  const TrackerTopology* tTopo) {
734  uint32_t detID = pixdet->geographicalId().rawId();
735  signal_map_type& theSignal = _signal[detID]; // check validity
736 
737  // Predefined efficiencies
738  float subdetEfficiency = 1.0;
739 
740  // setup the chip indices conversion
741  uint32_t Subid = DetId(detID).subdetId();
742  if (Subid == PixelSubdetector::PixelBarrel || Subid == StripSubdetector::TOB) { // barrel layers
743  uint32_t layerIndex = tTopo->pxbLayer(detID);
744  if (layerIndex - 1 < eff.barrel_efficiencies.size())
745  subdetEfficiency = eff.barrel_efficiencies[layerIndex - 1];
746  } else { // forward disks
747  uint32_t diskIndex = 2 * tTopo->pxfDisk(detID) - tTopo->pxfSide(detID);
748  if (diskIndex - 1 < eff.endcap_efficiencies.size())
749  subdetEfficiency = eff.endcap_efficiencies[diskIndex - 1];
750  }
751 
752  LogDebug("Phase2TrackerDigitizerAlgorithm") << " enter pixel_inefficiency " << subdetEfficiency;
753 
754  // Now loop again over pixels to kill some of them.
755  // Loop over hits, amplitude in electrons, channel = coded row,col
756  for (auto& s : theSignal) {
757  float rand = rengine_->flat();
758  if (rand > subdetEfficiency) {
759  // make amplitude =0
760  s.second.set(0.); // reset amplitude
761  }
762  }
763 }
764 
765 void Phase2TrackerDigitizerAlgorithm::initializeEvent(CLHEP::HepRandomEngine& eng) {
767  gaussDistribution_ = std::make_unique<CLHEP::RandGaussQ>(eng, 0., theNoiseInElectrons_);
768  }
769  // Threshold smearing with gaussian distribution:
770  if (addThresholdSmearing_) {
772  std::make_unique<CLHEP::RandGaussQ>(eng, theThresholdInE_Endcap_, theThresholdSmearing_Endcap_);
774  std::make_unique<CLHEP::RandGaussQ>(eng, theThresholdInE_Barrel_, theThresholdSmearing_Barrel_);
775  }
776  rengine_ = &eng;
777  _signal.clear();
778 }
779 
780 // =======================================================================================
781 //
782 // Set the drift direction accoring to the Bfield in local det-unit frame
783 // Works for both barrel and forward pixels.
784 // Replace the sign convention to fit M.Swartz's formulaes.
785 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
786 // parameter value
787 
789  const GlobalVector& bfield,
790  const DetId& detId) const {
791  Frame detFrame(pixdet->surface().position(), pixdet->surface().rotation());
792  LocalVector Bfield = detFrame.toLocal(bfield);
793 
794  float dir_x = 0.0;
795  float dir_y = 0.0;
796  float dir_z = 0.0;
797  float scale = 1.0;
798 
799  uint32_t detID = pixdet->geographicalId().rawId();
800  uint32_t Sub_detid = DetId(detID).subdetId();
801 
802  // Read Lorentz angle from DB:
803  if (use_LorentzAngle_DB_) {
804  bool isPixel = (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == PixelSubdetector::PixelEndcap);
805  if (isPixel) {
806  LogDebug("Phase2TrackerDigitizerAlgorithm") << "Read Lorentz angle from DB for Phase-2 IT";
807  } else {
808  LogDebug("Phase2TrackerDigitizerAlgorithm") << "Read Lorentz angle from DB for Phase-2 OT";
809  }
810 
811  float lorentzAngle =
813  float alpha2 = std::pow(lorentzAngle, 2);
814 
815  dir_x = -(lorentzAngle * Bfield.y() + alpha2 * Bfield.z() * Bfield.x());
816  dir_y = +(lorentzAngle * Bfield.x() - alpha2 * Bfield.z() * Bfield.y());
817  dir_z = -(1 + alpha2 * std::pow(Bfield.z(), 2));
818  scale = (1 + alpha2 * std::pow(Bfield.z(), 2));
819  } else {
820  // Read Lorentz angle from cfg file:
821  LogDebug("Phase2TrackerDigitizerAlgorithm") << "Read Lorentz angle from cfg file";
822  float alpha2_Endcap = 0.0;
823  float alpha2_Barrel = 0.0;
824  if (alpha2Order_) {
825  alpha2_Endcap = std::pow(tanLorentzAnglePerTesla_Endcap_, 2);
826  alpha2_Barrel = std::pow(tanLorentzAnglePerTesla_Barrel_, 2);
827  }
828 
829  if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) { // barrel layers
830  dir_x = -(tanLorentzAnglePerTesla_Barrel_ * Bfield.y() + alpha2_Barrel * Bfield.z() * Bfield.x());
831  dir_y = +(tanLorentzAnglePerTesla_Barrel_ * Bfield.x() - alpha2_Barrel * Bfield.z() * Bfield.y());
832  dir_z = -(1 + alpha2_Barrel * std::pow(Bfield.z(), 2));
833  scale = (1 + alpha2_Barrel * std::pow(Bfield.z(), 2));
834 
835  } else { // forward disks
836  dir_x = -(tanLorentzAnglePerTesla_Endcap_ * Bfield.y() + alpha2_Endcap * Bfield.z() * Bfield.x());
837  dir_y = +(tanLorentzAnglePerTesla_Endcap_ * Bfield.x() - alpha2_Endcap * Bfield.z() * Bfield.y());
838  dir_z = -(1 + alpha2_Endcap * std::pow(Bfield.z(), 2));
839  scale = (1 + alpha2_Endcap * std::pow(Bfield.z(), 2));
840  }
841  }
842 
843  LocalVector theDriftDirection = LocalVector(dir_x / scale, dir_y / scale, dir_z / scale);
844  LogDebug("Phase2TrackerDigitizerAlgorithm") << " The drift direction in local coordinate is " << theDriftDirection;
845  return theDriftDirection;
846 }
847 
848 // =============================================================================
849 
851  signal_map_type& theSignal = _signal[detID]; // check validity
852 
853  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
854  for (auto& s : theSignal) {
855  std::pair<int, int> ip;
856  if (pixelFlag_)
857  ip = PixelDigi::channelToPixel(s.first); //get pixel pos
858  else
859  ip = Phase2TrackerDigi::channelToPixel(s.first); //get pixel pos
860 
861  int row = ip.first; // X in row
862  int col = ip.second; // Y is in col
863  // transform to ROC index coordinates
864  if (theSiPixelGainCalibrationService_->isDead(detID, col, row))
865  s.second.set(0.); // reset amplitude
866  }
867 }
868 
869 // ==========================================================================
870 
872  bool isbad = false;
873  int detid = detID;
874  std::string Module;
875  for (auto const& det_m : deadModules_) {
876  int Dead_detID = det_m.getParameter<int>("Dead_detID");
877  Module = det_m.getParameter<std::string>("Module");
878  if (detid == Dead_detID) {
879  isbad = true;
880  break;
881  }
882  }
883 
884  if (!isbad)
885  return;
886 
887  signal_map_type& theSignal = _signal[detID]; // check validity
888  for (auto& s : theSignal) {
889  std::pair<int, int> ip;
890  if (pixelFlag_)
891  ip = PixelDigi::channelToPixel(s.first);
892  else
893  ip = Phase2TrackerDigi::channelToPixel(s.first); //get pixel pos
894 
895  if (Module == "whole" || (Module == "tbmA" && ip.first >= 80 && ip.first <= 159) ||
896  (Module == "tbmB" && ip.first <= 79))
897  s.second.set(0.);
898  }
899 }
900 // For premixing
901 void Phase2TrackerDigitizerAlgorithm::loadAccumulator(uint32_t detId, const std::map<int, float>& accumulator) {
902  auto& theSignal = _signal[detId];
903  // the input channel is always with PixelDigi definition
904  // if needed, that has to be converted to Phase2TrackerDigi convention
905  for (const auto& elem : accumulator) {
906  auto inserted = theSignal.emplace(elem.first, digitizerUtility::Ph2Amplitude(elem.second, nullptr));
907  if (!inserted.second) {
908  throw cms::Exception("LogicError") << "Signal was already set for DetId " << detId;
909  }
910  }
911 }
912 
914  std::map<int, digitizerUtility::DigiSimInfo>& digi_map,
915  const TrackerTopology* tTopo) {
916  uint32_t detID = pixdet->geographicalId().rawId();
917  auto it = _signal.find(detID);
918  if (it == _signal.end())
919  return;
920 
921  const signal_map_type& theSignal = _signal[detID];
922 
923  uint32_t Sub_detid = DetId(detID).subdetId();
924 
925  float theThresholdInE = 0.;
926  float theHIPThresholdInE = 0.;
927  // Define Threshold
928  if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) { // Barrel modules
929  theThresholdInE = addThresholdSmearing_ ? smearedThreshold_Barrel_->fire() // gaussian smearing
930  : theThresholdInE_Barrel_; // no smearing
931  theHIPThresholdInE = theHIPThresholdInE_Barrel_;
932  } else { // Forward disks modules
933  theThresholdInE = addThresholdSmearing_ ? smearedThreshold_Endcap_->fire() // gaussian smearing
934  : theThresholdInE_Endcap_; // no smearing
935  theHIPThresholdInE = theHIPThresholdInE_Endcap_;
936  }
937 
938  // if (addNoise) add_noise(pixdet, theThresholdInE/theNoiseInElectrons_); // generate noise
939  if (addNoise_)
940  add_noise(pixdet); // generate noise
941  if (addXtalk_)
942  add_cross_talk(pixdet);
943  if (addNoisyPixels_) {
944  float thresholdInNoiseUnits = 99.9;
946  thresholdInNoiseUnits = theThresholdInE / theNoiseInElectrons_;
947  add_noisy_cells(pixdet, thresholdInNoiseUnits);
948  }
949 
950  // Do only if needed
951  if (addPixelInefficiency_ && !theSignal.empty()) {
952  if (use_ineff_from_db_)
953  pixel_inefficiency_db(detID);
954  else
955  pixel_inefficiency(subdetEfficiencies_, pixdet, tTopo);
956  }
957  if (use_module_killing_) {
958  if (use_deadmodule_DB_) // remove dead modules using DB
959  module_killing_DB(pixdet);
960  else // remove dead modules using the list in cfg file
961  module_killing_conf(detID);
962  }
963 
964  // Digitize if the signal is greater than threshold
965  for (auto const& s : theSignal) {
966  const digitizerUtility::Ph2Amplitude& sig_data = s.second;
967  float signalInElectrons = sig_data.ampl();
968 
969  const auto& info_list = sig_data.simInfoList();
970  const digitizerUtility::SimHitInfo* hitInfo = nullptr;
971  if (!info_list.empty())
972  hitInfo = std::max_element(info_list.begin(), info_list.end())->second.get();
973 
974  if (isAboveThreshold(hitInfo, signalInElectrons, theThresholdInE)) { // check threshold
976  info.sig_tot = convertSignalToAdc(detID, signalInElectrons, theThresholdInE); // adc
977  info.ot_bit = signalInElectrons > theHIPThresholdInE ? true : false;
978  if (makeDigiSimLinks_) {
979  for (auto const& l : sig_data.simInfoList()) {
980  float charge_frac = l.first / signalInElectrons;
981  if (l.first > -5.0)
982  info.simInfoList.push_back({charge_frac, l.second.get()});
983  }
984  }
985  digi_map.insert({s.first, info});
986  }
987  }
988 }
989 //
990 // Scale the Signal using Dual Slope option
991 //
992 int Phase2TrackerDigitizerAlgorithm::convertSignalToAdc(uint32_t detID, float signal_in_elec, float threshold) {
993  int signal_in_adc;
994  int temp_signal;
995  const int max_limit = 10;
996  if (thePhase2ReadoutMode_ == 0) {
997  signal_in_adc = theAdcFullScale_;
998  } else {
999  if (thePhase2ReadoutMode_ == -1) {
1000  temp_signal = std::min(static_cast<int>(std::round(signal_in_elec / theElectronPerADC_)), theAdcFullScale_);
1001  } else {
1002  // calculate the kink point and the slope
1003  int dualslope_param = std::min(std::abs(thePhase2ReadoutMode_), max_limit);
1004  int kink_point = static_cast<int>(theAdcFullScale_ / 2) + 1;
1005  // C-ROC: first valid ToT code above threshold is 0b0000
1006  temp_signal = std::floor((signal_in_elec - threshold) / theElectronPerADC_);
1007  if (temp_signal > kink_point)
1008  temp_signal =
1009  static_cast<int>(std::floor((temp_signal - kink_point) / (pow(2, dualslope_param - 1)))) + kink_point;
1010  }
1011  signal_in_adc = std::min(temp_signal, theAdcFullScale_);
1012  LogTrace("Phase2TrackerDigitizerAlgorithm")
1013  << " DetId " << detID << " signal_in_elec " << signal_in_elec << " threshold " << threshold
1014  << " signal_above_thr " << signal_in_elec - threshold << " temp conversion "
1015  << std::floor((signal_in_elec - threshold) / theElectronPerADC_) + 1 << " signal after slope correction "
1016  << temp_signal << " signal_in_adc " << signal_in_adc;
1017  }
1018  return signal_in_adc;
1019 }
1020 
1022  constexpr float p1 = 12.5f;
1023  constexpr float p2 = 0.2733f;
1024  constexpr float p3 = 0.147f;
1025 
1026  auto xx = std::min(0.5f * x * x, p1);
1027  return 0.5f * (1.f - std::copysign(std::sqrt(1.f - unsafe_expf<4>(-xx * (1.f + p2 / (1.f + p3 * xx)))), x));
1028 }
virtual std::vector< digitizerUtility::EnergyDepositUnit > primary_ionization(const PSimHit &hit) const
virtual bool select_hit(const PSimHit &hit, double tCorr, double &sigScale) const
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned int pxbLayer(const DetId &id) const
static const TGPicture * info(bool iBackgroundIsBlack)
Local3DVector LocalVector
Definition: LocalVector.h:12
T perp() const
Definition: PV3DBase.h:69
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
void loadAccumulator(uint32_t detId, const std::map< int, float > &accumulator)
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:69
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
virtual int ncolumns() const =0
T z() const
Definition: PV3DBase.h:61
LocalVector driftDirection(const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const DetId &detId) const
virtual void initializeEvent(CLHEP::HepRandomEngine &eng)
virtual void induce_signal(std::vector< PSimHit >::const_iterator inputBegin, const PSimHit &hit, const size_t hitIndex, const size_t firstHitIndex, const uint32_t tofBin, const Phase2TrackerGeomDetUnit *pixdet, const std::vector< digitizerUtility::SignalPoint > &collection_points)
virtual int nrows() const =0
virtual void add_cross_talk(const Phase2TrackerGeomDetUnit *pixdet)
const SiPixelLorentzAngle * siPixelLorentzAngle_
T x() const
Definition: PV2DBase.h:43
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate_
const std::unique_ptr< SiPixelChargeReweightingAlgorithm > theSiPixelChargeReweightingAlgorithm_
std::map< int, digitizerUtility::Ph2Amplitude, std::less< int > > signal_map_type
virtual void add_noise(const Phase2TrackerGeomDetUnit *pixdet)
virtual void module_killing_conf(uint32_t detID)
#define LogTrace(id)
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:75
T y() const
Definition: PV2DBase.h:44
virtual float thickness() const =0
U second(std::pair< T, U > const &p)
const std::unique_ptr< SiPixelGainCalibrationOfflineSimService > theSiPixelGainCalibrationService_
virtual void add_noisy_cells(const Phase2TrackerGeomDetUnit *pixdet, float thePixelThreshold)
std::unique_ptr< CLHEP::RandGaussQ > gaussDistribution_
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
virtual void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, const size_t inputBeginGlobalIndex, const uint32_t tofBin, const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield)
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
virtual bool isAboveThreshold(const digitizerUtility::SimHitInfo *hitInfo, float charge, float thr) const
T sqrt(T t)
Definition: SSEVec.h:19
static PackedDigiType pixelToChannel(unsigned int row, unsigned int col)
Phase2TrackerDigitizerAlgorithm(const edm::ParameterSet &conf_common, const edm::ParameterSet &conf_specific, edm::ConsumesCollector iC)
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
constexpr float Bfield
Definition: Config.h:55
T mag() const
Definition: PV3DBase.h:64
unsigned int pxfDisk(const DetId &id) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int channel(const LocalPoint &p) const =0
double f[11][100]
virtual std::vector< float > fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments) const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
virtual void pixel_inefficiency_db(uint32_t detID)
static constexpr auto TOB
std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Endcap_
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
#define M_PI
std::vector< edm::ParameterSet > Parameters
constexpr double c_inv
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
Definition: DetId.h:17
unsigned int pxfSide(const DetId &id) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const PositionType & position() const
const SiPhase2OuterTrackerLorentzAngle * siPhase2OTLorentzAngle_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
virtual void digitize(const Phase2TrackerGeomDetUnit *pixdet, std::map< int, digitizerUtility::DigiSimInfo > &digi_map, const TrackerTopology *tTopo)
std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Barrel_
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi...
static std::pair< unsigned int, unsigned int > channelToPixel(unsigned int ch)
HLT enums.
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
col
Definition: cuy.py:1009
bool isPixel(HitType hitType)
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)
const RotationType & rotation() const
virtual void pixel_inefficiency(const SubdetEfficiencies &eff, const Phase2TrackerGeomDetUnit *pixdet, const TrackerTopology *tTopo)
virtual std::vector< digitizerUtility::SignalPoint > drift(const PSimHit &hit, const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const std::vector< digitizerUtility::EnergyDepositUnit > &ionization_points) const
virtual std::pair< float, float > pitch() const =0
Log< level::Warning, false > LogWarning
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
float getLorentzAngle(const uint32_t &) const
const std::vector< std::pair< float, std::unique_ptr< SimHitInfo > > > & simInfoList() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
virtual void module_killing_DB(const Phase2TrackerGeomDetUnit *pixdet)=0
int convertSignalToAdc(uint32_t detID, float signal_in_elec, float threshold)
#define LogDebug(id)
const Bounds & bounds() const
Definition: Surface.h:87
unsigned transform(const HcalDetId &id, unsigned transformCode)