CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Phase2TrackerDigitizerAlgorithm.cc
Go to the documentation of this file.
1 #include <typeinfo>
2 #include <iostream>
3 #include <cmath>
4 
6 
9 
10 #include "CLHEP/Random/RandGaussQ.h"
11 
16 
22 
31 
32 // Geometry
36 
37 using namespace edm;
38 using namespace sipixelobjects;
39 
40 namespace ph2tkdigialgo {
41  // Mass in MeV
42  constexpr double m_pion = 139.571;
43  constexpr double m_kaon = 493.677;
44  constexpr double m_electron = 0.511;
45  constexpr double m_muon = 105.658;
46  constexpr double m_proton = 938.272;
47 } // namespace ph2tkdigialgo
49  const edm::ParameterSet& conf_specific,
51  : _signal(),
52  makeDigiSimLinks_(conf_common.getUntrackedParameter<bool>("makeDigiSimLinks", true)),
53  use_ineff_from_db_(conf_specific.getParameter<bool>("Inefficiency_DB")),
54  use_module_killing_(conf_specific.getParameter<bool>("KillModules")), // boolean to kill or not modules
55  use_deadmodule_DB_(conf_specific.getParameter<bool>("DeadModules_DB")), // boolean to access dead modules from DB
56  // boolean to access Lorentz angle from DB
57  use_LorentzAngle_DB_(conf_specific.getParameter<bool>("LorentzAngle_DB")),
58 
59  // get dead module from cfg file
60  deadModules_(use_deadmodule_DB_ ? Parameters() : conf_specific.getParameter<Parameters>("DeadModules")),
61 
62  // Common pixel parameters
63  // These are parameters which are not likely to be changed
64  GeVperElectron_(3.61E-09), // 1 electron(3.61eV, 1keV(277e, mod 9/06 d.k.
65  alpha2Order_(conf_specific.getParameter<bool>("Alpha2Order")), // switch on/off of E.B effect
66  addXtalk_(conf_specific.getParameter<bool>("AddXTalk")),
67  // Interstrip Coupling - Not used in PixelDigitizerAlgorithm
68  interstripCoupling_(conf_specific.getParameter<double>("InterstripCoupling")),
69 
70  Sigma0_(conf_specific.getParameter<double>("SigmaZero")), // Charge diffusion constant 7->3.7
71  SigmaCoeff_(conf_specific.getParameter<double>("SigmaCoeff")), // delta in the diffusion across the strip pitch
72  // (set between 0 to 0.9, 0-->flat Sigma0, 1-->Sigma_min=0 & Sigma_max=2*Sigma0
73  // D.B.: Dist300 replaced by moduleThickness, may not work with partially depleted sensors but works otherwise
74  // Dist300(0.0300), // normalized to 300micron Silicon
75 
76  // Charge integration spread on the collection plane
77  clusterWidth_(conf_specific.getParameter<double>("ClusterWidth")),
78 
79  // Allowed modes of readout which has following values :
80  // 0 ---> Digital or binary readout
81  // -1 ---> Analog readout, current digitizer (Inner Pixel) (TDR version) with no threshold subtraction
82  // Analog readout with dual slope with the "second" slope being 1/2^(n-1) and threshold subtraction (n = 1, 2, 3,4)
83  thePhase2ReadoutMode_(conf_specific.getParameter<int>("Phase2ReadoutMode")),
84 
85  // ADC calibration 1adc count(135e.
86  // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
87  // Be careful, this parameter is also used in SiPixelDet.cc to
88  // calculate the noise in adc counts from noise in electrons.
89  // Both defaults should be the same.
90  theElectronPerADC_(conf_specific.getParameter<double>("ElectronPerAdc")),
91 
92  // ADC saturation value, 255(8bit adc.
93  theAdcFullScale_(conf_specific.getParameter<int>("AdcFullScale")),
94 
95  // Noise in electrons:
96  // Pixel cell noise, relevant for generating noisy pixels
97  theNoiseInElectrons_(conf_specific.getParameter<double>("NoiseInElectrons")),
98 
99  // Fill readout noise, including all readout chain, relevant for smearing
100  theReadoutNoise_(conf_specific.getParameter<double>("ReadoutNoiseInElec")),
101 
102  // Threshold in units of noise:
103  // thePixelThreshold(conf.getParameter<double>("ThresholdInNoiseUnits")),
104  // Pixel threshold in electron units.
105  theThresholdInE_Endcap_(conf_specific.getParameter<double>("ThresholdInElectrons_Endcap")),
106  theThresholdInE_Barrel_(conf_specific.getParameter<double>("ThresholdInElectrons_Barrel")),
107 
108  // Add threshold gaussian smearing:
109  theThresholdSmearing_Endcap_(conf_specific.getParameter<double>("ThresholdSmearing_Endcap")),
110  theThresholdSmearing_Barrel_(conf_specific.getParameter<double>("ThresholdSmearing_Barrel")),
111 
112  // Add HIP Threshold in electron units.
113  theHIPThresholdInE_Endcap_(conf_specific.getParameter<double>("HIPThresholdInElectrons_Endcap")),
114  theHIPThresholdInE_Barrel_(conf_specific.getParameter<double>("HIPThresholdInElectrons_Barrel")),
115 
116  // theTofCut 12.5, cut in particle TOD +/- 12.5ns
117  theTofLowerCut_(conf_specific.getParameter<double>("TofLowerCut")),
118  theTofUpperCut_(conf_specific.getParameter<double>("TofUpperCut")),
119 
120  // Get the Lorentz angle from the cfg file:
121  tanLorentzAnglePerTesla_Endcap_(
122  use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Endcap")),
123  tanLorentzAnglePerTesla_Barrel_(
124  use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Barrel")),
125 
126  // Add noise
127  addNoise_(conf_specific.getParameter<bool>("AddNoise")),
128 
129  // Add noisy pixels
130  addNoisyPixels_(conf_specific.getParameter<bool>("AddNoisyPixels")),
131 
132  // Fluctuate charge in track subsegments
133  fluctuateCharge_(conf_specific.getUntrackedParameter<bool>("FluctuateCharge", true)),
134 
135  // Control the pixel inefficiency
136  addPixelInefficiency_(conf_specific.getParameter<bool>("AddInefficiency")),
137 
138  // Add threshold gaussian smearing:
139  addThresholdSmearing_(conf_specific.getParameter<bool>("AddThresholdSmearing")),
140 
141  // Add some pseudo-red damage
142  pseudoRadDamage_(conf_specific.exists("PseudoRadDamage") ? conf_specific.getParameter<double>("PseudoRadDamage")
143  : double(0.0)),
144  pseudoRadDamageRadius_(conf_specific.exists("PseudoRadDamageRadius")
145  ? conf_specific.getParameter<double>("PseudoRadDamageRadius")
146  : double(0.0)),
147 
148  // delta cutoff in MeV, has to be same as in OSCAR(0.030/cmsim=1.0 MeV
149  // tMax(0.030), // In MeV.
150  // tMax(conf.getUntrackedParameter<double>("DeltaProductionCut",0.030)),
151  tMax_(conf_common.getParameter<double>("DeltaProductionCut")),
152 
153  badPixels_(conf_specific.getParameter<Parameters>("CellsToKill")),
154 
155  fluctuate_(fluctuateCharge_ ? std::make_unique<SiG4UniversalFluctuation>() : nullptr),
156  theNoiser_(addNoise_ ? std::make_unique<GaussianTailNoiseGenerator>() : nullptr),
157  theSiPixelGainCalibrationService_(
158  use_ineff_from_db_ ? std::make_unique<SiPixelGainCalibrationOfflineSimService>(conf_specific, iC) : nullptr),
159  subdetEfficiencies_(conf_specific) {
160  LogDebug("Phase2TrackerDigitizerAlgorithm")
161  << "Phase2TrackerDigitizerAlgorithm constructed\n"
162  << "Configuration parameters:\n"
163  << "Threshold/Gain = "
164  << "threshold in electron Endcap = " << theThresholdInE_Endcap_
165  << "\nthreshold in electron Barrel = " << theThresholdInE_Barrel_ << " ElectronPerADC " << theElectronPerADC_
166  << " ADC Scale (in bits) " << theAdcFullScale_ << " The delta cut-off is set to " << tMax_ << " pix-inefficiency "
168 }
169 
171  LogDebug("Phase2TrackerDigitizerAlgorithm") << "Phase2TrackerDigitizerAlgorithm deleted";
172 }
173 
175  barrel_efficiencies = conf.getParameter<std::vector<double> >("EfficiencyFactors_Barrel");
176  endcap_efficiencies = conf.getParameter<std::vector<double> >("EfficiencyFactors_Endcap");
177 }
178 void Phase2TrackerDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
179  std::vector<PSimHit>::const_iterator inputEnd,
180  const size_t inputBeginGlobalIndex,
181  const uint32_t tofBin,
182  const Phase2TrackerGeomDetUnit* pixdet,
183  const GlobalVector& bfield) {
184  // produce SignalPoint's for all SimHit's in detector
185  // Loop over hits
186  uint32_t detId = pixdet->geographicalId().rawId();
187  size_t simHitGlobalIndex = inputBeginGlobalIndex; // This needs to be stored to create the digi-sim link later
188 
189  // find the relevant hits
190  std::vector<PSimHit> matchedSimHits;
191  std::copy_if(inputBegin, inputEnd, std::back_inserter(matchedSimHits), [detId](auto const& hit) -> bool {
192  return hit.detUnitId() == detId;
193  });
194  // loop over a much reduced set of SimHits
195  for (auto const& hit : matchedSimHits) {
196  LogDebug("Phase2DigitizerAlgorithm") << hit.particleType() << " " << hit.pabs() << " " << hit.energyLoss() << " "
197  << hit.tof() << " " << hit.trackId() << " " << hit.processType() << " "
198  << hit.detUnitId() << hit.entryPoint() << " " << hit.exitPoint();
199  double signalScale = 1.0;
200  // fill collection_points for this SimHit, indpendent of topology
201  if (select_hit(hit, (pixdet->surface().toGlobal(hit.localPosition()).mag() * c_inv), signalScale)) {
202  const auto& ionization_points = primary_ionization(hit); // fills ionization_points
203 
204  // transforms ionization_points -> collection_points
205  const auto& collection_points = drift(hit, pixdet, bfield, ionization_points);
206 
207  // compute induced signal on readout elements and add to _signal
208  // hit needed only for SimHit<-->Digi link
209  induce_signal(hit, simHitGlobalIndex, tofBin, pixdet, collection_points);
210  }
211  ++simHitGlobalIndex;
212  }
213 }
214 // =================================================================
215 //
216 // Generate primary ionization along the track segment.
217 // Divide the track into small sub-segments
218 //
219 // =================================================================
220 std::vector<DigitizerUtility::EnergyDepositUnit> Phase2TrackerDigitizerAlgorithm::primary_ionization(
221  const PSimHit& hit) const {
222  // Straight line approximation for trajectory inside active media
223  constexpr float SegmentLength = 0.0010; // in cm (10 microns)
224  // Get the 3D segment direction vector
225  LocalVector direction = hit.exitPoint() - hit.entryPoint();
226 
227  float eLoss = hit.energyLoss(); // Eloss in GeV
228  float length = direction.mag(); // Track length in Silicon
229 
230  int NumberOfSegments = static_cast<int>(length / SegmentLength); // Number of segments
231  if (NumberOfSegments < 1)
232  NumberOfSegments = 1;
233  LogDebug("Phase2TrackerDigitizerAlgorithm")
234  << "enter primary_ionzation " << NumberOfSegments << " shift = " << hit.exitPoint().x() - hit.entryPoint().x()
235  << " " << hit.exitPoint().y() - hit.entryPoint().y() << " " << hit.exitPoint().z() - hit.entryPoint().z() << " "
236  << hit.particleType() << " " << hit.pabs();
237 
238  std::vector<float> elossVector;
239  elossVector.reserve(NumberOfSegments);
240  if (fluctuateCharge_) {
241  // Generate fluctuated charge points
242  elossVector = fluctuateEloss(hit.particleType(), hit.pabs(), eLoss, length, NumberOfSegments);
243  } else {
244  float averageEloss = eLoss / NumberOfSegments;
245  elossVector.resize(NumberOfSegments, averageEloss);
246  }
247 
248  std::vector<DigitizerUtility::EnergyDepositUnit> ionization_points;
249  ionization_points.reserve(NumberOfSegments); // set size
250  // loop over segments
251  for (size_t i = 0; i < elossVector.size(); ++i) {
252  // Divide the segment into equal length subsegments
253  Local3DPoint point = hit.entryPoint() + ((i + 0.5) / NumberOfSegments) * direction;
254  float energy = elossVector[i] / GeVperElectron_; // Convert charge to elec.
255 
256  DigitizerUtility::EnergyDepositUnit edu(energy, point); // define position,energy point
257  ionization_points.push_back(edu); // save
258  LogDebug("Phase2TrackerDigitizerAlgorithm")
259  << i << " " << edu.x() << " " << edu.y() << " " << edu.z() << " " << edu.energy();
260  }
261  return ionization_points;
262 }
263 //==============================================================================
264 //
265 // Fluctuate the charge comming from a small (10um) track segment.
266 // Use the G4 routine. For mip pions for the moment.
267 //
268 //==============================================================================
270  int pid, float particleMomentum, float eloss, float length, int NumberOfSegs) const {
271  double particleMass = ph2tkdigialgo::m_pion; // Mass in MeV, assume pion
272  pid = std::abs(pid);
273  if (pid != 211) { // Mass in MeV
274  if (pid == 11)
275  particleMass = ph2tkdigialgo::m_electron;
276  else if (pid == 13)
277  particleMass = ph2tkdigialgo::m_muon;
278  else if (pid == 321)
279  particleMass = ph2tkdigialgo::m_kaon;
280  else if (pid == 2212)
281  particleMass = ph2tkdigialgo::m_proton;
282  }
283  // What is the track segment length.
284  float segmentLength = length / NumberOfSegs;
285 
286  // Generate charge fluctuations.
287  float sum = 0.;
288  double segmentEloss = (1000. * eloss) / NumberOfSegs; //eloss in MeV
289  std::vector<float> elossVector;
290  elossVector.reserve(NumberOfSegs);
291  for (int i = 0; i < NumberOfSegs; ++i) {
292  // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
293  // track segment length in mm, segment eloss in MeV
294  // Returns fluctuated eloss in MeV
295  double deltaCutoff = tMax_; // the cutoff is sometimes redefined inside, so fix it.
296  float de = fluctuate_->SampleFluctuations(particleMomentum * 1000.,
297  particleMass,
298  deltaCutoff,
299  segmentLength * 10.,
300  segmentEloss,
301  rengine_) /
302  1000.; //convert to GeV
303  elossVector.push_back(de);
304  sum += de;
305  }
306  if (sum > 0.) { // if fluctuations give eloss>0.
307  // Rescale to the same total eloss
308  float ratio = eloss / sum;
310  std::begin(elossVector), std::end(elossVector), std::begin(elossVector), [&ratio](auto const& c) -> float {
311  return c * ratio;
312  }); // use a simple lambda expression
313  } else { // if fluctuations gives 0 eloss
314  float averageEloss = eloss / NumberOfSegs;
315  elossVector.resize(NumberOfSegs, averageEloss);
316  }
317  return elossVector;
318 }
319 
320 // ======================================================================
321 //
322 // Drift the charge segments to the sensor surface (collection plane)
323 // Include the effect of E-field and B-field
324 //
325 // =====================================================================
326 std::vector<DigitizerUtility::SignalPoint> Phase2TrackerDigitizerAlgorithm::drift(
327  const PSimHit& hit,
328  const Phase2TrackerGeomDetUnit* pixdet,
329  const GlobalVector& bfield,
330  const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points) const {
331  LogDebug("Phase2TrackerDigitizerAlgorithm") << "enter drift ";
332 
333  std::vector<DigitizerUtility::SignalPoint> collection_points;
334  collection_points.reserve(ionization_points.size()); // set size
335  LocalVector driftDir = DriftDirection(pixdet, bfield, hit.detUnitId()); // get the charge drift direction
336  if (driftDir.z() == 0.) {
337  LogWarning("Phase2TrackerDigitizerAlgorithm") << " pxlx: drift in z is zero ";
338  return collection_points;
339  }
340 
341  float TanLorenzAngleX = driftDir.x(); // tangent of Lorentz angle
342  float TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
343  float dir_z = driftDir.z(); // The z drift direction
344  float CosLorenzAngleX = 1. / std::sqrt(1. + std::pow(TanLorenzAngleX, 2)); // cosine to estimate the path length
345  float CosLorenzAngleY = 1.;
346  if (alpha2Order_) {
347  TanLorenzAngleY = driftDir.y();
348  CosLorenzAngleY = 1. / std::sqrt(1. + std::pow(TanLorenzAngleY, 2)); // cosine
349  }
350 
351  float moduleThickness = pixdet->specificSurface().bounds().thickness();
352  float stripPitch = pixdet->specificTopology().pitch().first;
353 
354  LogDebug("Phase2TrackerDigitizerAlgorithm")
355  << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY << " " << CosLorenzAngleX << " "
356  << CosLorenzAngleY << " " << moduleThickness * TanLorenzAngleX << " " << driftDir;
357 
358  for (auto const& val : ionization_points) {
359  // position
360  float SegX = val.x(), SegY = val.y(), SegZ = val.z();
361 
362  // Distance from the collection plane
363  // DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
364  // Include explixitely the E drift direction (for CMS dir_z=-1)
365 
366  // Distance between charge generation and collection
367  float driftDistance = moduleThickness / 2. - (dir_z * SegZ); // Drift to -z
368 
369  if (driftDistance < 0.)
370  driftDistance = 0.;
371  else if (driftDistance > moduleThickness)
372  driftDistance = moduleThickness;
373 
374  // Assume full depletion now, partial depletion will come later.
375  float XDriftDueToMagField = driftDistance * TanLorenzAngleX;
376  float YDriftDueToMagField = driftDistance * TanLorenzAngleY;
377 
378  // Shift cloud center
379  float CloudCenterX = SegX + XDriftDueToMagField;
380  float CloudCenterY = SegY + YDriftDueToMagField;
381 
382  // Calculate how long is the charge drift path
383  // Actual Drift Lentgh
384  float driftLength =
385  std::sqrt(std::pow(driftDistance, 2) + std::pow(XDriftDueToMagField, 2) + std::pow(YDriftDueToMagField, 2));
386 
387  // What is the charge diffusion after this path
388  // Sigma0=0.00037 is for 300um thickness (make sure moduleThickness is in [cm])
389  float Sigma = std::sqrt(driftLength / moduleThickness) * Sigma0_ * moduleThickness / 0.0300;
390  // D.B.: sigmaCoeff=0 means no modulation
391  if (SigmaCoeff_)
392  Sigma *= (SigmaCoeff_ * std::pow(cos(SegX * M_PI / stripPitch), 2) + 1);
393  // NB: divided by 4 to get a periodicity of stripPitch
394 
395  // Project the diffusion sigma on the collection plane
396  float Sigma_x = Sigma / CosLorenzAngleX;
397  float Sigma_y = Sigma / CosLorenzAngleY;
398 
399  // Insert a charge loss due to Rad Damage here
400  float energyOnCollector = val.energy(); // The energy that reaches the collector
401 
402  // pseudoRadDamage
403  if (pseudoRadDamage_) {
404  float moduleRadius = pixdet->surface().position().perp();
405  if (moduleRadius <= pseudoRadDamageRadius_) {
406  float kValue = pseudoRadDamage_ / std::pow(moduleRadius, 2);
407  energyOnCollector *= exp(-1 * kValue * driftDistance / moduleThickness);
408  }
409  }
410  LogDebug("Phase2TrackerDigitizerAlgorithm")
411  << "Dift DistanceZ = " << driftDistance << " module thickness = " << moduleThickness
412  << " Start Energy = " << val.energy() << " Energy after loss= " << energyOnCollector;
413  DigitizerUtility::SignalPoint sp(CloudCenterX, CloudCenterY, Sigma_x, Sigma_y, hit.tof(), energyOnCollector);
414 
415  // Load the Charge distribution parameters
416  collection_points.push_back(sp);
417  }
418  return collection_points;
419 }
420 
421 // ====================================================================
422 //
423 // Induce the signal on the collection plane of the active sensor area.
425  const PSimHit& hit,
426  const size_t hitIndex,
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, " << topol->pitch().first << " " << topol->pitch().second;
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  float corr_time = hit.tof() - pixdet->surface().toGlobal(hit.localPosition()).mag() * c_inv;
571  for (auto const& hit_s : hit_signal) {
572  int chan = hit_s.first;
573  theSignal[chan] +=
574  (makeDigiSimLinks_ ? DigitizerUtility::Amplitude(hit_s.second, &hit, hit_s.second, corr_time, hitIndex, tofBin)
575  : DigitizerUtility::Amplitude(hit_s.second, nullptr, hit_s.second));
576  }
577 }
578 // ======================================================================
579 //
580 // Add electronic noise to pixel charge
581 //
582 // ======================================================================
584  uint32_t detID = pixdet->geographicalId().rawId();
585  signal_map_type& theSignal = _signal[detID];
586  for (auto& s : theSignal) {
587  float noise = gaussDistribution_->fire();
588  if ((s.second.ampl() + noise) < 0.)
589  s.second.set(0);
590  else
591  s.second += noise;
592  }
593 }
594 // ======================================================================
595 //
596 // Add Cross-talk contribution
597 //
598 // ======================================================================
600  uint32_t detID = pixdet->geographicalId().rawId();
601  signal_map_type& theSignal = _signal[detID];
602  signal_map_type signalNew;
603  const Phase2TrackerTopology* topol = &pixdet->specificTopology();
604  int numRows = topol->nrows();
605 
606  for (auto& s : theSignal) {
607  float signalInElectrons = s.second.ampl(); // signal in electrons
608 
609  std::pair<int, int> hitChan;
610  if (pixelFlag_)
611  hitChan = PixelDigi::channelToPixel(s.first);
612  else
613  hitChan = Phase2TrackerDigi::channelToPixel(s.first);
614 
615  float signalInElectrons_Xtalk = signalInElectrons * interstripCoupling_;
616  // subtract the charge which will be shared
617  s.second.set(signalInElectrons - signalInElectrons_Xtalk);
618 
619  if (hitChan.first != 0) {
620  auto XtalkPrev = std::make_pair(hitChan.first - 1, hitChan.second);
621  int chanXtalkPrev = pixelFlag_ ? PixelDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second)
622  : Phase2TrackerDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second);
623  signalNew.emplace(chanXtalkPrev, DigitizerUtility::Amplitude(signalInElectrons_Xtalk, nullptr, -1.0));
624  }
625  if (hitChan.first < numRows - 1) {
626  auto XtalkNext = std::make_pair(hitChan.first + 1, hitChan.second);
627  int chanXtalkNext = pixelFlag_ ? PixelDigi::pixelToChannel(XtalkNext.first, XtalkNext.second)
628  : Phase2TrackerDigi::pixelToChannel(XtalkNext.first, XtalkNext.second);
629  signalNew.emplace(chanXtalkNext, DigitizerUtility::Amplitude(signalInElectrons_Xtalk, nullptr, -1.0));
630  }
631  }
632  for (auto const& l : signalNew) {
633  int chan = l.first;
634  auto iter = theSignal.find(chan);
635  if (iter != theSignal.end()) {
636  theSignal[chan] += l.second.ampl();
637  } else {
638  theSignal.emplace(chan, DigitizerUtility::Amplitude(l.second.ampl(), nullptr, -1.0));
639  }
640  }
641 }
642 
643 // ======================================================================
644 //
645 // Add noise on non-hit cells
646 //
647 // ======================================================================
649  uint32_t detID = pixdet->geographicalId().rawId();
650  signal_map_type& theSignal = _signal[detID];
651  const Phase2TrackerTopology* topol = &pixdet->specificTopology();
652 
653  int numColumns = topol->ncolumns(); // det module number of cols&rows
654  int numRows = topol->nrows();
655 
656  int numberOfPixels = numRows * numColumns;
657  std::map<int, float, std::less<int> > otherPixels;
658 
659  theNoiser_->generate(numberOfPixels,
660  thePixelThreshold, //thr. in un. of nois
661  theNoiseInElectrons_, // noise in elec.
662  otherPixels,
663  rengine_);
664 
665  LogDebug("Phase2TrackerDigitizerAlgorithm")
666  << " Add noisy pixels " << numRows << " " << numColumns << " " << theNoiseInElectrons_ << " "
667  << theThresholdInE_Endcap_ << " " << theThresholdInE_Barrel_ << " " << numberOfPixels << " "
668  << otherPixels.size();
669 
670  // Add noisy pixels
671  for (auto const& el : otherPixels) {
672  int iy = el.first / numRows;
673  if (iy < 0 || iy > numColumns - 1)
674  LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in iy " << iy;
675 
676  int ix = el.first - iy * numRows;
677  if (ix < 0 || ix > numRows - 1)
678  LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in ix " << ix;
679 
681 
682  LogDebug("Phase2TrackerDigitizerAlgorithm")
683  << " Storing noise = " << el.first << " " << el.second << " " << ix << " " << iy << " " << chan;
684 
685  if (theSignal[chan] == 0)
686  theSignal[chan] = DigitizerUtility::Amplitude(el.second, nullptr, -1.);
687  }
688 }
689 // ============================================================================
690 //
691 // Simulate the readout inefficiencies.
692 // Delete a selected number of single pixels, dcols and rocs.
694  const Phase2TrackerGeomDetUnit* pixdet,
695  const TrackerTopology* tTopo) {
696  uint32_t detID = pixdet->geographicalId().rawId();
697  signal_map_type& theSignal = _signal[detID]; // check validity
698 
699  // Predefined efficiencies
700  float subdetEfficiency = 1.0;
701 
702  // setup the chip indices conversion
703  uint32_t Subid = DetId(detID).subdetId();
704  if (Subid == PixelSubdetector::PixelBarrel || Subid == StripSubdetector::TOB) { // barrel layers
705  uint32_t layerIndex = tTopo->pxbLayer(detID);
706  if (layerIndex - 1 < eff.barrel_efficiencies.size())
707  subdetEfficiency = eff.barrel_efficiencies[layerIndex - 1];
708  } else { // forward disks
709  uint32_t diskIndex = 2 * tTopo->pxfDisk(detID) - tTopo->pxfSide(detID);
710  if (diskIndex - 1 < eff.endcap_efficiencies.size())
711  subdetEfficiency = eff.endcap_efficiencies[diskIndex - 1];
712  }
713 
714  LogDebug("Phase2TrackerDigitizerAlgorithm") << " enter pixel_inefficiency " << subdetEfficiency;
715 
716  // Now loop again over pixels to kill some of them.
717  // Loop over hits, amplitude in electrons, channel = coded row,col
718  for (auto& s : theSignal) {
719  float rand = rengine_->flat();
720  if (rand > subdetEfficiency) {
721  // make amplitude =0
722  s.second.set(0.); // reset amplitude
723  }
724  }
725 }
726 void Phase2TrackerDigitizerAlgorithm::initializeEvent(CLHEP::HepRandomEngine& eng) {
728  gaussDistribution_ = std::make_unique<CLHEP::RandGaussQ>(eng, 0., theReadoutNoise_);
729  }
730  // Threshold smearing with gaussian distribution:
731  if (addThresholdSmearing_) {
733  std::make_unique<CLHEP::RandGaussQ>(eng, theThresholdInE_Endcap_, theThresholdSmearing_Endcap_);
735  std::make_unique<CLHEP::RandGaussQ>(eng, theThresholdInE_Barrel_, theThresholdSmearing_Barrel_);
736  }
737  rengine_ = &eng;
738  _signal.clear();
739 }
740 
741 // =======================================================================================
742 //
743 // Set the drift direction accoring to the Bfield in local det-unit frame
744 // Works for both barrel and forward pixels.
745 // Replace the sign convention to fit M.Swartz's formulaes.
746 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
747 // parameter value
748 
750  const GlobalVector& bfield,
751  const DetId& detId) const {
752  Frame detFrame(pixdet->surface().position(), pixdet->surface().rotation());
753  LocalVector Bfield = detFrame.toLocal(bfield);
754 
755  float dir_x = 0.0;
756  float dir_y = 0.0;
757  float dir_z = 0.0;
758  float scale = 1.0;
759 
760  uint32_t detID = pixdet->geographicalId().rawId();
761  uint32_t Sub_detid = DetId(detID).subdetId();
762 
763  // Read Lorentz angle from DB:
764  if (use_LorentzAngle_DB_) {
765  bool isPixel = (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == PixelSubdetector::PixelEndcap);
766 
767  float lorentzAngle =
769  float alpha2 = std::pow(lorentzAngle, 2);
770 
771  dir_x = -(lorentzAngle * Bfield.y() + alpha2 * Bfield.z() * Bfield.x());
772  dir_y = +(lorentzAngle * Bfield.x() - alpha2 * Bfield.z() * Bfield.y());
773  dir_z = -(1 + alpha2 * std::pow(Bfield.z(), 2));
774  scale = (1 + alpha2 * std::pow(Bfield.z(), 2));
775  } else {
776  // Read Lorentz angle from cfg file:
777  float alpha2_Endcap = 0.0;
778  float alpha2_Barrel = 0.0;
779  if (alpha2Order_) {
780  alpha2_Endcap = std::pow(tanLorentzAnglePerTesla_Endcap_, 2);
781  alpha2_Barrel = std::pow(tanLorentzAnglePerTesla_Barrel_, 2);
782  }
783 
784  if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) { // barrel layers
785  dir_x = -(tanLorentzAnglePerTesla_Barrel_ * Bfield.y() + alpha2_Barrel * Bfield.z() * Bfield.x());
786  dir_y = +(tanLorentzAnglePerTesla_Barrel_ * Bfield.x() - alpha2_Barrel * Bfield.z() * Bfield.y());
787  dir_z = -(1 + alpha2_Barrel * std::pow(Bfield.z(), 2));
788  scale = (1 + alpha2_Barrel * std::pow(Bfield.z(), 2));
789 
790  } else { // forward disks
791  dir_x = -(tanLorentzAnglePerTesla_Endcap_ * Bfield.y() + alpha2_Endcap * Bfield.z() * Bfield.x());
792  dir_y = +(tanLorentzAnglePerTesla_Endcap_ * Bfield.x() - alpha2_Endcap * Bfield.z() * Bfield.y());
793  dir_z = -(1 + alpha2_Endcap * std::pow(Bfield.z(), 2));
794  scale = (1 + alpha2_Endcap * std::pow(Bfield.z(), 2));
795  }
796  }
797 
798  LocalVector theDriftDirection = LocalVector(dir_x / scale, dir_y / scale, dir_z / scale);
799  LogDebug("Phase2TrackerDigitizerAlgorithm") << " The drift direction in local coordinate is " << theDriftDirection;
800  return theDriftDirection;
801 }
802 
803 // =============================================================================
804 
806  signal_map_type& theSignal = _signal[detID]; // check validity
807 
808  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
809  for (auto& s : theSignal) {
810  std::pair<int, int> ip;
811  if (pixelFlag_)
812  ip = PixelDigi::channelToPixel(s.first); //get pixel pos
813  else
814  ip = Phase2TrackerDigi::channelToPixel(s.first); //get pixel pos
815 
816  int row = ip.first; // X in row
817  int col = ip.second; // Y is in col
818  // transform to ROC index coordinates
819  if (theSiPixelGainCalibrationService_->isDead(detID, col, row))
820  s.second.set(0.); // reset amplitude
821  }
822 }
823 
824 // ==========================================================================
825 
827  bool isbad = false;
828  int detid = detID;
829  std::string Module;
830  for (auto const& det_m : deadModules_) {
831  int Dead_detID = det_m.getParameter<int>("Dead_detID");
832  Module = det_m.getParameter<std::string>("Module");
833  if (detid == Dead_detID) {
834  isbad = true;
835  break;
836  }
837  }
838 
839  if (!isbad)
840  return;
841 
842  signal_map_type& theSignal = _signal[detID]; // check validity
843  for (auto& s : theSignal) {
844  std::pair<int, int> ip;
845  if (pixelFlag_)
846  ip = PixelDigi::channelToPixel(s.first);
847  else
848  ip = Phase2TrackerDigi::channelToPixel(s.first); //get pixel pos
849 
850  if (Module == "whole" || (Module == "tbmA" && ip.first >= 80 && ip.first <= 159) ||
851  (Module == "tbmB" && ip.first <= 79))
852  s.second.set(0.);
853  }
854 }
855 // ==========================================================================
857  bool isbad = false;
858  uint32_t detID = pixdet->geographicalId().rawId();
859  int ncol = pixdet->specificTopology().ncolumns();
860  if (ncol < 0)
861  return;
862  std::vector<SiPixelQuality::disabledModuleType> disabledModules = siPixelBadModule_->getBadComponentList();
863 
865  for (size_t id = 0; id < disabledModules.size(); id++) {
866  if (detID == disabledModules[id].DetID) {
867  isbad = true;
868  badmodule = disabledModules[id];
869  break;
870  }
871  }
872 
873  if (!isbad)
874  return;
875 
876  signal_map_type& theSignal = _signal[detID]; // check validity
877  if (badmodule.errorType == 0) { // this is a whole dead module.
878  for (auto& s : theSignal)
879  s.second.set(0.); // reset amplitude
880  } else { // all other module types: half-modules and single ROCs.
881  // Get Bad ROC position:
882  // follow the example of getBadRocPositions in CondFormats/SiPixelObjects/src/SiPixelQuality.cc
883  std::vector<GlobalPixel> badrocpositions;
884  for (size_t j = 0; j < static_cast<size_t>(ncol); j++) {
885  if (siPixelBadModule_->IsRocBad(detID, j)) {
886  std::vector<CablingPathToDetUnit> path = fedCablingMap_->pathToDetUnit(detID);
887  for (auto const& p : path) {
888  const PixelROC* myroc = fedCablingMap_->findItem(p);
889  if (myroc->idInDetUnit() == j) {
890  LocalPixel::RocRowCol local = {39, 25}; //corresponding to center of ROC row, col
891  GlobalPixel global = myroc->toGlobal(LocalPixel(local));
892  badrocpositions.push_back(global);
893  break;
894  }
895  }
896  }
897  }
898 
899  for (auto& s : theSignal) {
900  std::pair<int, int> ip;
901  if (pixelFlag_)
902  ip = PixelDigi::channelToPixel(s.first);
903  else
905 
906  for (auto const& p : badrocpositions) {
907  for (auto& k : badPixels_) {
908  if (p.row == k.getParameter<int>("row") && ip.first == k.getParameter<int>("row") &&
909  std::abs(ip.second - p.col) < k.getParameter<int>("col")) {
910  s.second.set(0.);
911  }
912  }
913  }
914  }
915  }
916 }
917 
918 // For premixing
919 void Phase2TrackerDigitizerAlgorithm::loadAccumulator(uint32_t detId, const std::map<int, float>& accumulator) {
920  auto& theSignal = _signal[detId];
921  // the input channel is always with PixelDigi definition
922  // if needed, that has to be converted to Phase2TrackerDigi convention
923  for (const auto& elem : accumulator) {
924  auto inserted = theSignal.emplace(elem.first, DigitizerUtility::Amplitude(elem.second, nullptr));
925  if (!inserted.second) {
926  throw cms::Exception("LogicError") << "Signal was already set for DetId " << detId;
927  }
928  }
929 }
930 
932  std::map<int, DigitizerUtility::DigiSimInfo>& digi_map,
933  const TrackerTopology* tTopo) {
934  uint32_t detID = pixdet->geographicalId().rawId();
935  auto it = _signal.find(detID);
936  if (it == _signal.end())
937  return;
938 
939  const signal_map_type& theSignal = _signal[detID];
940 
941  uint32_t Sub_detid = DetId(detID).subdetId();
942 
943  float theThresholdInE = 0.;
944  float theHIPThresholdInE = 0.;
945  // Define Threshold
946  if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) { // Barrel modules
947  theThresholdInE = addThresholdSmearing_ ? smearedThreshold_Barrel_->fire() // gaussian smearing
948  : theThresholdInE_Barrel_; // no smearing
949  theHIPThresholdInE = theHIPThresholdInE_Barrel_;
950  } else { // Forward disks modules
951  theThresholdInE = addThresholdSmearing_ ? smearedThreshold_Endcap_->fire() // gaussian smearing
952  : theThresholdInE_Endcap_; // no smearing
953  theHIPThresholdInE = theHIPThresholdInE_Endcap_;
954  }
955 
956  // if (addNoise) add_noise(pixdet, theThresholdInE/theNoiseInElectrons_); // generate noise
957  if (addNoise_)
958  add_noise(pixdet); // generate noise
959  if (addXtalk_)
960  add_cross_talk(pixdet);
961  if (addNoisyPixels_)
962  add_noisy_cells(pixdet, theHIPThresholdInE / theElectronPerADC_);
963 
964  // Do only if needed
965  if (addPixelInefficiency_ && !theSignal.empty()) {
966  if (use_ineff_from_db_)
967  pixel_inefficiency_db(detID);
968  else
969  pixel_inefficiency(subdetEfficiencies_, pixdet, tTopo);
970  }
971  if (use_module_killing_) {
972  if (use_deadmodule_DB_) // remove dead modules using DB
973  module_killing_DB(pixdet);
974  else // remove dead modules using the list in cfg file
975  module_killing_conf(detID);
976  }
977 
978  // Digitize if the signal is greater than threshold
979  for (auto const& s : theSignal) {
980  const DigitizerUtility::Amplitude& sig_data = s.second;
981  float signalInElectrons = sig_data.ampl();
982 
983  const auto& info_list = sig_data.simInfoList();
984  const DigitizerUtility::SimHitInfo* hitInfo = nullptr;
985  if (!info_list.empty())
986  hitInfo = std::max_element(info_list.begin(), info_list.end())->second.get();
987 
988  if (isAboveThreshold(hitInfo, signalInElectrons, theThresholdInE)) { // check threshold
990  info.sig_tot = convertSignalToAdc(detID, signalInElectrons, theThresholdInE); // adc
991  info.ot_bit = signalInElectrons > theHIPThresholdInE ? true : false;
992  if (makeDigiSimLinks_) {
993  for (auto const& l : sig_data.simInfoList()) {
994  float charge_frac = l.first / signalInElectrons;
995  if (l.first > -5.0)
996  info.simInfoList.push_back({charge_frac, l.second.get()});
997  }
998  }
999  digi_map.insert({s.first, info});
1000  }
1001  }
1002 }
1003 //
1004 // Scale the Signal using Dual Slope option
1005 //
1006 int Phase2TrackerDigitizerAlgorithm::convertSignalToAdc(uint32_t detID, float signal_in_elec, float threshold) {
1007  int signal_in_adc;
1008  int temp_signal;
1009  const int max_limit = 10;
1010  if (thePhase2ReadoutMode_ == 0) {
1011  signal_in_adc = theAdcFullScale_;
1012  } else {
1013  if (thePhase2ReadoutMode_ == -1) {
1014  temp_signal = std::min(static_cast<int>(signal_in_elec / theElectronPerADC_), theAdcFullScale_);
1015  } else {
1016  // calculate the kink point and the slope
1017  int dualslope_param = std::min(std::abs(thePhase2ReadoutMode_), max_limit);
1018  int kink_point = static_cast<int>(theAdcFullScale_ / 2) + 1;
1019  // C-ROC: first valid ToT code above threshold is 0b0000
1020  temp_signal = std::floor((signal_in_elec - threshold) / theElectronPerADC_);
1021  if (temp_signal > kink_point)
1022  temp_signal = std::floor((temp_signal - kink_point) / (pow(2, dualslope_param - 1))) + kink_point;
1023  }
1024  signal_in_adc = std::min(temp_signal, theAdcFullScale_);
1025  LogTrace("Phase2TrackerDigitizerAlgorithm")
1026  << " DetId " << detID << " signal_in_elec " << signal_in_elec << " threshold " << threshold
1027  << " signal_above_thr " << signal_in_elec - threshold << " temp conversion "
1028  << std::floor((signal_in_elec - threshold) / theElectronPerADC_) + 1 << " signal after slope correction "
1029  << temp_signal << " signal_in_adc " << signal_in_adc;
1030  }
1031  return signal_in_adc;
1032 }
1034  constexpr float p1 = 12.5f;
1035  constexpr float p2 = 0.2733f;
1036  constexpr float p3 = 0.147f;
1037 
1038  auto xx = std::min(0.5f * x * x, p1);
1039  return 0.5f * (1.f - std::copysign(std::sqrt(1.f - unsafe_expf<4>(-xx * (1.f + p2 / (1.f + p3 * xx)))), x));
1040 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
virtual void module_killing_DB(const Phase2TrackerGeomDetUnit *pixdet)
LocalVector DriftDirection(const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const DetId &detId) const
static const TGPicture * info(bool iBackgroundIsBlack)
Local3DVector LocalVector
Definition: LocalVector.h:12
const edm::EventSetup & c
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:76
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
void loadAccumulator(uint32_t detId, const std::map< int, float > &accumulator)
T y() const
Definition: PV2DBase.h:44
virtual bool select_hit(const PSimHit &hit, double tCorr, double &sigScale) const
T perp() const
Definition: PV3DBase.h:69
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:69
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
uint16_t *__restrict__ id
virtual int ncolumns() const =0
const TString p2
Definition: fwPaths.cc:13
virtual bool isAboveThreshold(const DigitizerUtility::SimHitInfo *hitInfo, float charge, float thr) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
unsigned int pxfDisk(const DetId &id) const
tuple chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, &quot;NDC&quot;) lumi.SetBorderSize( 0 ) lumi...
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T y() const
Definition: PV3DBase.h:60
virtual void initializeEvent(CLHEP::HepRandomEngine &eng)
virtual int nrows() const =0
const Bounds & bounds() const
Definition: Surface.h:87
std::vector< sipixelobjects::CablingPathToDetUnit > pathToDetUnit(uint32_t rawDetId) const final
virtual void add_cross_talk(const Phase2TrackerGeomDetUnit *pixdet)
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
const SiPixelLorentzAngle * siPixelLorentzAngle_
bool IsRocBad(const uint32_t &detid, const short &rocNb) const
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate_
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
virtual void add_noise(const Phase2TrackerGeomDetUnit *pixdet)
virtual void module_killing_conf(uint32_t detID)
identify pixel inside single ROC
Definition: LocalPixel.h:7
#define LogTrace(id)
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:75
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
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_
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.
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
T mag() const
Definition: PV3DBase.h:64
Local3DPoint localPosition() const
Definition: PSimHit.h:52
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_
T z() const
Definition: PV3DBase.h:61
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:37
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
const std::vector< std::pair< float, std::unique_ptr< SimHitInfo > > > & simInfoList() const
const sipixelobjects::PixelROC * findItem(const sipixelobjects::CablingPathToDetUnit &path) const final
const TString p1
Definition: fwPaths.cc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int channel(const LocalPoint &p) const =0
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
virtual void induce_signal(const PSimHit &hit, const size_t hitIndex, const uint32_t tofBin, const Phase2TrackerGeomDetUnit *pixdet, const std::vector< DigitizerUtility::SignalPoint > &collection_points)
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
virtual void pixel_inefficiency_db(uint32_t detID)
const std::vector< disabledModuleType > getBadComponentList() const
T min(T a, T b)
Definition: MathUtil.h:58
static constexpr auto TOB
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Endcap_
virtual void digitize(const Phase2TrackerGeomDetUnit *pixdet, std::map< int, DigitizerUtility::DigiSimInfo > &digi_map, const TrackerTopology *tTopo)
std::map< int, DigitizerUtility::Amplitude, std::less< int > > signal_map_type
#define M_PI
std::vector< edm::ParameterSet > Parameters
std::vector< std::pair< float, SimHitInfo * > > simInfoList
constexpr double c_inv
unsigned int pxbLayer(const DetId &id) const
__shared__ int noise
Definition: DetId.h:17
const SiPhase2OuterTrackerLorentzAngle * siPhase2OTLorentzAngle_
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
row and collumn in ROC representation
Definition: LocalPixel.h:13
std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Barrel_
virtual std::vector< DigitizerUtility::EnergyDepositUnit > primary_ionization(const PSimHit &hit) const
static std::pair< unsigned int, unsigned int > channelToPixel(unsigned int ch)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
int particleType() const
Definition: PSimHit.h:89
string end
Definition: dataset.py:937
bool isPixel(HitType hitType)
float getLorentzAngle(const uint32_t &) const
unsigned int pxfSide(const DetId &id) const
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
virtual std::vector< float > fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments) const
int col
Definition: cuy.py:1009
T x() const
Definition: PV2DBase.h:43
T x() const
Definition: PV3DBase.h:59
const PositionType & position() const
float getLorentzAngle(const uint32_t &) const
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
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
unsigned int detUnitId() const
Definition: PSimHit.h:97
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
int convertSignalToAdc(uint32_t detID, float signal_in_elec, float threshold)
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:55
#define LogDebug(id)
unsigned transform(const HcalDetId &id, unsigned transformCode)