CMS 3D CMS Logo

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