CMS 3D CMS Logo

Phase2TrackerDigitizerAlgorithm.cc
Go to the documentation of this file.
1 #include<typeinfo>
2 #include <iostream>
3 #include <cmath>
4 
6 
10 
12 #include "CLHEP/Random/RandGaussQ.h"
13 #include "CLHEP/Random/RandFlat.h"
14 
15 //#include "PixelIndices.h"
20 
26 
41 
48 
49 // Geometry
54 
56 using namespace edm;
57 using namespace sipixelobjects;
58 
60  const edm::ParameterSet& conf_specific,
61  CLHEP::HepRandomEngine& eng):
62  _signal(),
63  makeDigiSimLinks_(conf_specific.getUntrackedParameter<bool>("makeDigiSimLinks", true)),
64  use_ineff_from_db_(conf_specific.getParameter<bool>("Inefficiency_DB")),
65  use_module_killing_(conf_specific.getParameter<bool>("KillModules")), // boolean to kill or not modules
66  use_deadmodule_DB_(conf_specific.getParameter<bool>("DeadModules_DB")), // boolean to access dead modules from DB
67  use_LorentzAngle_DB_(conf_specific.getParameter<bool>("LorentzAngle_DB")), // boolean to access Lorentz angle from DB
68 
69  DeadModules(use_deadmodule_DB_ ? Parameters() : conf_specific.getParameter<Parameters>("DeadModules")), // get dead module from cfg file
70 
71  // Common pixel parameters
72  // These are parameters which are not likely to be changed
73  GeVperElectron(3.61E-09), // 1 electron(3.61eV, 1keV(277e, mod 9/06 d.k.
74  alpha2Order(conf_specific.getParameter<bool>("Alpha2Order")), // switch on/off of E.B effect
75  addXtalk(conf_specific.getParameter<bool>("AddXTalk")),
76  interstripCoupling(conf_specific.getParameter<double>("InterstripCoupling")), // Interstrip Coupling
77 
78  Sigma0(conf_specific.getParameter<double>("SigmaZero")), // Charge diffusion constant 7->3.7
79  SigmaCoeff(conf_specific.getParameter<double>("SigmaCoeff")), // delta in the diffusion across the strip pitch
80  // (set between 0 to 0.9, 0-->flat Sigma0, 1-->Sigma_min=0 & Sigma_max=2*Sigma0
81  // D.B.: Dist300 replaced by moduleThickness, may not work with partially depleted sensors but works otherwise
82  // Dist300(0.0300), // normalized to 300micron Silicon
83 
84  ClusterWidth(conf_specific.getParameter<double>("ClusterWidth")), // Charge integration spread on the collection plane
85 
86  doDigitalReadout(conf_specific.getParameter<bool>("DigitalReadout")), // Flag to decide analog or digital readout
87 
88  // ADC calibration 1adc count(135e.
89  // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
90  // Be careful, this parameter is also used in SiPixelDet.cc to
91  // calculate the noise in adc counts from noise in electrons.
92  // Both defaults should be the same.
93  theElectronPerADC(conf_specific.getParameter<double>("ElectronPerAdc")),
94 
95  // ADC saturation value, 255(8bit adc.
96  theAdcFullScale(conf_specific.getParameter<int>("AdcFullScale")),
97 
98  // Noise in electrons:
99  // Pixel cell noise, relevant for generating noisy pixels
100  theNoiseInElectrons(conf_specific.getParameter<double>("NoiseInElectrons")),
101 
102  // Fill readout noise, including all readout chain, relevant for smearing
103  theReadoutNoise(conf_specific.getParameter<double>("ReadoutNoiseInElec")),
104 
105  // Threshold in units of noise:
106  // thePixelThreshold(conf.getParameter<double>("ThresholdInNoiseUnits")),
107  // Pixel threshold in electron units.
108  theThresholdInE_Endcap(conf_specific.getParameter<double>("ThresholdInElectrons_Endcap")),
109  theThresholdInE_Barrel(conf_specific.getParameter<double>("ThresholdInElectrons_Barrel")),
110 
111  // Add threshold gaussian smearing:
112  theThresholdSmearing_Endcap(conf_specific.getParameter<double>("ThresholdSmearing_Endcap")),
113  theThresholdSmearing_Barrel(conf_specific.getParameter<double>("ThresholdSmearing_Barrel")),
114 
115  // Add HIP Threshold in electron units.
116  theHIPThresholdInE_Endcap(conf_specific.getParameter<double>("HIPThresholdInElectrons_Endcap")),
117  theHIPThresholdInE_Barrel(conf_specific.getParameter<double>("HIPThresholdInElectrons_Barrel")),
118 
119  // theTofCut 12.5, cut in particle TOD +/- 12.5ns
120  theTofLowerCut(conf_specific.getParameter<double>("TofLowerCut")),
121  theTofUpperCut(conf_specific.getParameter<double>("TofUpperCut")),
122 
123  // Get the Lorentz angle from the cfg file:
124  tanLorentzAnglePerTesla_Endcap(use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Endcap")),
125  tanLorentzAnglePerTesla_Barrel(use_LorentzAngle_DB_ ? 0.0 : conf_specific.getParameter<double>("TanLorentzAnglePerTesla_Barrel")),
126 
127  // Add noise
128  addNoise(conf_specific.getParameter<bool>("AddNoise")),
129 
130  // Add noisy pixels
131  addNoisyPixels(conf_specific.getParameter<bool>("AddNoisyPixels")),
132 
133  // Fluctuate charge in track subsegments
134  fluctuateCharge(conf_specific.getUntrackedParameter<bool>("FluctuateCharge",true)),
135 
136  // Control the pixel inefficiency
137  AddPixelInefficiency(conf_specific.getParameter<bool>("AddInefficiency")),
138 
139  // Add threshold gaussian smearing:
140  addThresholdSmearing(conf_specific.getParameter<bool>("AddThresholdSmearing")),
141 
142  // Add some pseudo-red damage
143  pseudoRadDamage(conf_specific.exists("PseudoRadDamage")?conf_specific.getParameter<double>("PseudoRadDamage"):double(0.0)),
144  pseudoRadDamageRadius(conf_specific.exists("PseudoRadDamageRadius")?conf_specific.getParameter<double>("PseudoRadDamageRadius"):double(0.0)),
145 
146  // delta cutoff in MeV, has to be same as in OSCAR(0.030/cmsim=1.0 MeV
147  // tMax(0.030), // In MeV.
148  // tMax(conf.getUntrackedParameter<double>("DeltaProductionCut",0.030)),
149  tMax(conf_common.getParameter<double>("DeltaProductionCut")),
150 
151  badPixels(conf_specific.getParameter<std::vector<edm::ParameterSet> >("CellsToKill")),
152  fluctuate(fluctuateCharge ? new SiG4UniversalFluctuation() : nullptr),
153  theNoiser(addNoise ? new GaussianTailNoiseGenerator() : nullptr),
154  theSiPixelGainCalibrationService_(use_ineff_from_db_ ? new SiPixelGainCalibrationOfflineSimService(conf_specific) : nullptr),
155  subdetEfficiencies_(conf_specific),
156  flatDistribution_((addNoise || AddPixelInefficiency || fluctuateCharge || addThresholdSmearing) ? new CLHEP::RandFlat(eng, 0., 1.) : nullptr),
157  gaussDistribution_((addNoise || AddPixelInefficiency || fluctuateCharge || addThresholdSmearing) ? new CLHEP::RandGaussQ(eng, 0., theReadoutNoise) : nullptr),
158  // Threshold smearing with gaussian distribution:
159  smearedThreshold_Endcap_(addThresholdSmearing ? new CLHEP::RandGaussQ(eng, theThresholdInE_Endcap , theThresholdSmearing_Endcap) : nullptr),
160  smearedThreshold_Barrel_(addThresholdSmearing ? new CLHEP::RandGaussQ(eng, theThresholdInE_Barrel , theThresholdSmearing_Barrel) : nullptr),
161  rengine_(&eng)
162 {
163  LogInfo("Phase2TrackerDigitizerAlgorithm") << "Phase2TrackerDigitizerAlgorithm constructed\n"
164  << "Configuration parameters:\n"
165  << "Threshold/Gain = "
166  << "threshold in electron Endcap = "
168  << "\nthreshold in electron Barrel = "
170  << " " << theElectronPerADC << " " << theAdcFullScale
171  << " The delta cut-off is set to " << tMax
172  << " pix-inefficiency " << AddPixelInefficiency;
173 }
174 
176  LogDebug("Phase2TrackerDigitizerAlgorithm") << "Phase2TrackerDigitizerAlgorithm deleted";
177 }
178 
180  barrel_efficiencies = conf.getParameter< std::vector<double> >("EfficiencyFactors_Barrel");
181  endcap_efficiencies = conf.getParameter< std::vector<double> >("EfficiencyFactors_Endcap");
182 }
183 // =================================================================
184 //
185 // Generate primary ionization along the track segment.
186 // Divide the track into small sub-segments
187 //
188 // =================================================================
190  std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points) const {
191  // Straight line approximation for trajectory inside active media
192  const float SegmentLength = 0.0010; // in cm (10 microns)
193  float energy;
194 
195  // Get the 3D segment direction vector
196  LocalVector direction = hit.exitPoint() - hit.entryPoint();
197 
198  float eLoss = hit.energyLoss(); // Eloss in GeV
199  float length = direction.mag(); // Track length in Silicon
200 
201  int NumberOfSegments = int (length / SegmentLength); // Number of segments
202  if (NumberOfSegments < 1) NumberOfSegments = 1;
203  LogDebug("Phase2TrackerDigitizerAlgorithm")
204  << "enter primary_ionzation " << NumberOfSegments
205  << " shift = "
206  << hit.exitPoint().x() - hit.entryPoint().x() << " "
207  << hit.exitPoint().y() - hit.entryPoint().y() << " "
208  << hit.exitPoint().z() - hit.entryPoint().z() << " "
209  << hit.particleType()
210  << " " << hit.pabs();
211 
212  std::vector<float> elossVector; // Eloss vector
213  elossVector.reserve(NumberOfSegments);
214  if (fluctuateCharge) {
215  int pid = hit.particleType();
216  // int pid=211; // assume it is a pion
217 
218  float momentum = hit.pabs();
219  // Generate fluctuated charge points
220  fluctuateEloss(pid, momentum, eLoss, length, NumberOfSegments, elossVector);
221  }
222  ionization_points.reserve(NumberOfSegments); // set size
223 
224  // loop over segments
225  for (int i = 0; i != NumberOfSegments; ++i) {
226  // Divide the segment into equal length subsegments
227  Local3DPoint point = hit.entryPoint() + float((i+0.5)/NumberOfSegments) * direction;
228  if (fluctuateCharge)
229  energy = elossVector[i]/GeVperElectron; // Convert charge to elec.
230  else
231  energy = hit.energyLoss()/GeVperElectron/float(NumberOfSegments);
232 
233  DigitizerUtility::EnergyDepositUnit edu(energy, point); // define position,energy point
234  ionization_points.push_back(edu); // save
235  LogDebug("Phase2TrackerDigitizerAlgorithm")
236  << i << " " << ionization_points[i].x() << " "
237  << ionization_points[i].y() << " "
238  << ionization_points[i].z() << " "
239  << ionization_points[i].energy();
240  }
241 }
242 //==============================================================================
243 //
244 // Fluctuate the charge comming from a small (10um) track segment.
245 // Use the G4 routine. For mip pions for the moment.
246 //
247 //==============================================================================
249  float particleMomentum,
250  float eloss,
251  float length,
252  int NumberOfSegs,
253  std::vector<float> & elossVector) const {
254 
255  // Get dedx for this track
256  //float dedx;
257  //if( length > 0.) dedx = eloss/length;
258  //else dedx = eloss;
259 
260  double particleMass = 139.6; // Mass in MeV, Assume pion
261  pid = std::abs(pid);
262  if (pid != 211) { // Mass in MeV
263  if (pid == 11) particleMass = 0.511;
264  else if (pid == 13) particleMass = 105.7;
265  else if (pid == 321) particleMass = 493.7;
266  else if (pid == 2212) particleMass = 938.3;
267  }
268  // What is the track segment length.
269  float segmentLength = length/NumberOfSegs;
270 
271  // Generate charge fluctuations.
272  float de = 0.;
273  float sum = 0.;
274  double segmentEloss = (1000. * eloss)/NumberOfSegs; //eloss in MeV
275  for (int i = 0; i < NumberOfSegs; ++i) {
276  // material,*, momentum,energy,*, *, mass
277  //myglandz_(14.,segmentLength,2.,2.,dedx,de,0.14);
278  // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
279  // track segment length in mm, segment eloss in MeV
280  // Returns fluctuated eloss in MeV
281  double deltaCutoff = tMax; // the cutoff is sometimes redefined inside, so fix it.
282  de = fluctuate->SampleFluctuations(static_cast<double>(particleMomentum*1000.),
283  particleMass, deltaCutoff,
284  static_cast<double>(segmentLength*10.),
285  segmentEloss,
286  rengine_ )/1000.; //convert to GeV
287  elossVector.push_back(de);
288  sum += de;
289  }
290  if (sum > 0.) { // If fluctuations give eloss>0.
291  // Rescale to the same total eloss
292  float ratio = eloss/sum;
293  for (int ii = 0; ii < NumberOfSegs; ++ii) elossVector[ii] = ratio*elossVector[ii];
294  }
295  else { // If fluctuations gives 0 eloss
296  float averageEloss = eloss/NumberOfSegs;
297  for (int ii = 0; ii < NumberOfSegs; ++ii) elossVector[ii] = averageEloss;
298  }
299 
300 }
301 
302 // ======================================================================
303 //
304 // Drift the charge segments to the sensor surface (collection plane)
305 // Include the effect of E-field and B-field
306 //
307 // =====================================================================
309  const Phase2TrackerGeomDetUnit* pixdet,
310  const GlobalVector& bfield,
311  const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points,
312  std::vector<DigitizerUtility::SignalPoint>& collection_points) const {
313  LogDebug("Phase2TrackerDigitizerAlgorithm") << "enter drift ";
314 
315  collection_points.resize(ionization_points.size()); // set size
316  LocalVector driftDir = DriftDirection(pixdet, bfield, hit.detUnitId()); // get the charge drift direction
317  if (driftDir.z() == 0.) {
318  LogWarning("Phase2TrackerDigitizerAlgorithm") << " pxlx: drift in z is zero ";
319  return;
320  }
321 
322  float TanLorenzAngleX,
323  TanLorenzAngleY,
324  dir_z,
325  CosLorenzAngleX,
326  CosLorenzAngleY;
327  if (alpha2Order) {
328  TanLorenzAngleX = driftDir.x(); // tangen of Lorentz angle
329  TanLorenzAngleY = driftDir.y();
330  dir_z = driftDir.z(); // The z drift direction
331  CosLorenzAngleX = 1./std::sqrt(1. + TanLorenzAngleX * TanLorenzAngleX); // cosine
332  CosLorenzAngleY = 1./std::sqrt(1. + TanLorenzAngleY * TanLorenzAngleY); // cosine;
333  }
334  else {
335  TanLorenzAngleX = driftDir.x();
336  TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
337  dir_z = driftDir.z(); // The z drift direction
338  CosLorenzAngleX = 1./std::sqrt(1. + TanLorenzAngleX * TanLorenzAngleX); // cosine to estimate the path length
339  CosLorenzAngleY = 1.;
340  }
341 
342  float moduleThickness = pixdet->specificSurface().bounds().thickness();
343  float stripPitch = pixdet->specificTopology().pitch().first;
344 
345  LogDebug("Phase2TrackerDigitizerAlgorithm")
346  << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY <<" "
347  << CosLorenzAngleX << " " << CosLorenzAngleY << " "
348  << moduleThickness*TanLorenzAngleX << " " << driftDir;
349 
350  float Sigma_x = 1.; // Charge spread
351  float Sigma_y = 1.;
352  float DriftDistance; // Distance between charge generation and collection
353  float DriftLength; // Actual Drift Lentgh
354  float Sigma;
355 
356  for (unsigned int i = 0; i != ionization_points.size(); ++i) {
357  float SegX, SegY, SegZ; // position
358  SegX = ionization_points[i].x();
359  SegY = ionization_points[i].y();
360  SegZ = ionization_points[i].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  DriftDistance = moduleThickness/2. - (dir_z * SegZ); // Drift to -z
366 
367  if (DriftDistance < 0.)
368  DriftDistance = 0.;
369  else if (DriftDistance > moduleThickness)
370  DriftDistance = moduleThickness;
371 
372  // Assume full depletion now, partial depletion will come later.
373  float XDriftDueToMagField = DriftDistance * TanLorenzAngleX;
374  float YDriftDueToMagField = DriftDistance * TanLorenzAngleY;
375 
376  // Shift cloud center
377  float CloudCenterX = SegX + XDriftDueToMagField;
378  float CloudCenterY = SegY + YDriftDueToMagField;
379 
380  // Calculate how long is the charge drift path
381  DriftLength = std::sqrt(DriftDistance*DriftDistance +
382  XDriftDueToMagField*XDriftDueToMagField +
383  YDriftDueToMagField*YDriftDueToMagField);
384 
385  // What is the charge diffusion after this path
386  // Sigma0=0.00037 is for 300um thickness (make sure moduleThickness is in [cm])
387  Sigma = std::sqrt(DriftLength/moduleThickness) * (Sigma0 * moduleThickness/0.0300);
388  // D.B.: sigmaCoeff=0 means no modulation
389  if (SigmaCoeff) Sigma *= (SigmaCoeff * cos(SegX*M_PI/stripPitch) * cos(SegX*M_PI/stripPitch) + 1);
390  // NB: divided by 4 to get a periodicity of stripPitch
391 
392  // Project the diffusion sigma on the collection plane
393  Sigma_x = Sigma / CosLorenzAngleX;
394  Sigma_y = Sigma / CosLorenzAngleY;
395 
396  // Insert a charge loss due to Rad Damage here
397  float energyOnCollector = ionization_points[i].energy(); // The energy that reaches the collector
398 
399  // pseudoRadDamage
400  if (pseudoRadDamage >= 0.001) {
401  float moduleRadius = pixdet->surface().position().perp();
402  if (moduleRadius <= pseudoRadDamageRadius) {
403  float kValue = pseudoRadDamage/(moduleRadius * moduleRadius);
404  energyOnCollector = energyOnCollector * exp(-1 * kValue * DriftDistance/moduleThickness);
405  }
406  }
407  LogDebug("Phase2TrackerDigitizerAlgorithm")
408  << "Dift DistanceZ = " << DriftDistance << " module thickness = " << moduleThickness
409  << " Start Energy = " << ionization_points[i].energy() << " Energy after loss= " << energyOnCollector;
410  DigitizerUtility::SignalPoint sp(CloudCenterX, CloudCenterY, Sigma_x, Sigma_y, hit.tof(), energyOnCollector);
411  // Load the Charge distribution parameters
412  collection_points[i] = sp;
413  }
414 }
415 
416 // ====================================================================
417 //
418 // Induce the signal on the collection plane of the active sensor area.
420  const size_t hitIndex,
421  const unsigned int tofBin,
422  const Phase2TrackerGeomDetUnit* pixdet,
423  const std::vector<DigitizerUtility::SignalPoint>& collection_points) {
424 
425  // X - Rows, Left-Right, 160, (1.6cm) for barrel
426  // Y - Columns, Down-Up, 416, (6.4cm)
427  const Phase2TrackerTopology* topol = &pixdet->specificTopology();
428  uint32_t detID = pixdet->geographicalId().rawId();
429  signal_map_type& theSignal = _signal[detID];
430 
431  LogDebug("Phase2TrackerDigitizerAlgorithm")
432  << " enter induce_signal, "
433  << topol->pitch().first << " " << topol->pitch().second; //OK
434 
435  // local map to store pixels hit by 1 Hit.
436  using hit_map_type = std::map<int, float, std::less<int> >;
437  hit_map_type hit_signal;
438 
439  // map to store pixel integrals in the x and in the y directions
440  std::map<int, float, std::less<int> > x,y;
441 
442  // Assign signals to readout channels and store sorted by channel number
443  int iseg = 0;
444  float ESum = 0.0;
445 
446  // Iterate over collection points on the collection plane
447  for (auto const & v : collection_points) {
448  iseg++;
449  float CloudCenterX = v.position().x(); // Charge position in x
450  float CloudCenterY = v.position().y(); // in y
451  float SigmaX = v.sigma_x(); // Charge spread in x
452  float SigmaY = v.sigma_y(); // in y
453  float Charge = v.amplitude(); // Charge amplitude
454 
455  LogDebug("Phase2TrackerDigitizerAlgorithm")
456  << " cloud " << v.position().x() << " " << v.position().y() << " "
457  << v.sigma_x() << " " << v.sigma_y() << " " << v.amplitude();
458 
459  // Find the maximum cloud spread in 2D plane , assume 3*sigma
460  float CloudRight = CloudCenterX + ClusterWidth * SigmaX;
461  float CloudLeft = CloudCenterX - ClusterWidth * SigmaX;
462  float CloudUp = CloudCenterY + ClusterWidth * SigmaY;
463  float CloudDown = CloudCenterY - ClusterWidth * SigmaY;
464 
465  // Define 2D cloud limit points
466  LocalPoint PointRightUp = LocalPoint(CloudRight,CloudUp);
467  LocalPoint PointLeftDown = LocalPoint(CloudLeft,CloudDown);
468 
469  // This points can be located outside the sensor area.
470  // The conversion to measurement point does not check for that
471  // so the returned pixel index might be wrong (outside range).
472  // We rely on the limits check below to fix this.
473  // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
474  // PIXEL ARE IS LOST, it should not be collected.
475 
476  // Convert the 2D points to pixel indices
477  MeasurementPoint mp = topol->measurementPosition(PointRightUp ); //OK
478 
479  int IPixRightUpX = int(floor( mp.x()));
480  int IPixRightUpY = int(floor( mp.y()));
481 
482  LogDebug("Phase2TrackerDigitizerAlgorithm") << " right-up " << PointRightUp << " "
483  << mp.x() << " " << mp.y() << " "
484  << IPixRightUpX << " " << IPixRightUpY ;
485 
486  mp = topol->measurementPosition(PointLeftDown); // OK
487 
488  int IPixLeftDownX = int(floor( mp.x()));
489  int IPixLeftDownY = int(floor( mp.y()));
490 
491  LogDebug("Phase2TrackerDigitizerAlgorithm") << " left-down " << PointLeftDown << " "
492  << 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  x.clear(); // clear temporary integration array
505  y.clear();
506 
507  // First integrate cahrge strips in x
508  int ix; // TT for compatibility
509  for (ix = IPixLeftDownX; ix <= IPixRightUpX; ix++) { // loop over x index
510  float xUB, xLB, UpperBound, LowerBound;
511 
512  // Why is set to 0 if ix=0, does it meen that we accept charge
513  // outside the sensor?
514  if (ix == 0 || SigmaX==0.) // skip for surface segemnts
515  LowerBound = 0.;
516  else {
517  mp = MeasurementPoint( float(ix), 0.0);
518  xLB = topol->localPosition(mp).x();
519  LowerBound = 1-calcQ((xLB-CloudCenterX)/SigmaX);
520  }
521 
522  if (ix == numRows-1 || SigmaX == 0.)
523  UpperBound = 1.;
524  else {
525  mp = MeasurementPoint(float(ix+1), 0.0);
526  xUB = topol->localPosition(mp).x();
527  UpperBound = 1. - calcQ((xUB-CloudCenterX)/SigmaX);
528  }
529  float TotalIntegrationRange = UpperBound - LowerBound; // get strip
530  x[ix] = TotalIntegrationRange; // save strip integral
531  }
532 
533  // Now integarte strips in y
534  int iy; // TT for compatibility
535  for (iy = IPixLeftDownY; iy <= IPixRightUpY; iy++) { //loope over y ind
536  float yUB, yLB, UpperBound, LowerBound;
537 
538  if (iy == 0 || SigmaY==0.)
539  LowerBound = 0.;
540  else {
541  mp = MeasurementPoint(0.0, float(iy));
542  yLB = topol->localPosition(mp).y();
543  LowerBound = 1. - calcQ((yLB-CloudCenterY)/SigmaY);
544  }
545 
546  if (iy == numColumns-1 || SigmaY==0. )
547  UpperBound = 1.;
548  else {
549  mp = MeasurementPoint(0.0, float(iy+1));
550  yUB = topol->localPosition(mp).y();
551  UpperBound = 1. - calcQ((yUB-CloudCenterY)/SigmaY);
552  }
553 
554  float TotalIntegrationRange = UpperBound - LowerBound;
555  y[iy] = TotalIntegrationRange; // save strip integral
556  }
557 
558  // Get the 2D charge integrals by folding x and y strips
559  int chan;
560  for (ix = IPixLeftDownX; ix <= IPixRightUpX; ix++) { // loop over x index
561  for (iy = IPixLeftDownY; iy <= IPixRightUpY; iy++) { //loope over y ind
562  float ChargeFraction = Charge*x[ix]*y[iy];
563  if (ChargeFraction > 0.) {
564  chan = (pixelFlag) ? PixelDigi::pixelToChannel(ix, iy)
565  : Phase2TrackerDigi::pixelToChannel(ix, iy); // Get index
566  // Load the amplitude
567  hit_signal[chan] += ChargeFraction;
568  }
569 
570  mp = MeasurementPoint(float(ix), float(iy));
571  LocalPoint lp = topol->localPosition(mp);
572  chan = topol->channel(lp);
573 
574  LogDebug("Phase2TrackerDigitizerAlgorithm")
575  << " pixel " << ix << " " << iy << " - "<<" "
576  << chan << " " << ChargeFraction<<" "
577  << mp.x() << " " << mp.y() <<" "
578  << lp.x() << " " << lp.y() << " " // givex edge position
579  << chan; // edge belongs to previous ?
580  ESum += ChargeFraction;
581  }
582  }
583  }
584  // Fill the global map with all hit pixels from this event
585  for (auto const & hit_s : hit_signal) {
586  int chan = hit_s.first;
587  theSignal[chan] += (makeDigiSimLinks_ ? DigitizerUtility::Amplitude( hit_s.second, &hit, hit_s.second, hitIndex, tofBin)
588  : DigitizerUtility::Amplitude( hit_s.second, nullptr, hit_s.second) ) ;
589  }
590 }
591 // ======================================================================
592 //
593 // Add electronic noise to pixel charge
594 //
595 // ======================================================================
597  float thePixelThreshold) {
598  LogDebug("Phase2TrackerDigitizerAlgorithm") << " enter add_noise " << theNoiseInElectrons;
599  uint32_t detID = pixdet->geographicalId().rawId();
600  signal_map_type& theSignal = _signal[detID];
601  const Phase2TrackerTopology* topol = &pixdet->specificTopology();
602  int numColumns = topol->ncolumns(); // det module number of cols&rows
603  int numRows = topol->nrows();
604 
605  // First add Readout noise to hit cells
606  for (auto & s : theSignal) {
607  float noise = gaussDistribution_->fire();
608  if ((s.second.ampl() + noise) < 0.)
609  s.second.set(0);
610  else
611  s.second += noise;
612  }
613 
614  if (addXtalk) {
615  signal_map_type signalNew;
616  for (auto const & s : theSignal) {
617  float signalInElectrons = s.second.ampl(); // signal in electrons
618  std::pair<int,int> hitChan;
619  if (pixelFlag) hitChan = PixelDigi::channelToPixel(s.first);
620  else hitChan = Phase2TrackerDigi::channelToPixel(s.first);
621 
622  float signalInElectrons_Xtalk = signalInElectrons * interstripCoupling;
623 
624  if (hitChan.first != 0) {
625  std::pair<int,int> XtalkPrev = std::pair<int,int>(hitChan.first-1, hitChan.second);
626  int chanXtalkPrev = (pixelFlag) ? PixelDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second)
627  : Phase2TrackerDigi::pixelToChannel(XtalkPrev.first, XtalkPrev.second);
628  signalNew.insert(std::pair<int,DigitizerUtility::Amplitude>(chanXtalkPrev, DigitizerUtility::Amplitude(signalInElectrons_Xtalk, nullptr, -1.0)));
629  }
630  if (hitChan.first < (numRows-1)) {
631  std::pair<int,int> XtalkNext = std::pair<int,int>(hitChan.first+1, hitChan.second);
632  int chanXtalkNext = (pixelFlag) ? PixelDigi::pixelToChannel(XtalkNext.first, XtalkNext.second)
633  : Phase2TrackerDigi::pixelToChannel(XtalkNext.first, XtalkNext.second);
634  signalNew.insert(std::pair<int,DigitizerUtility::Amplitude>(chanXtalkNext, DigitizerUtility::Amplitude(signalInElectrons_Xtalk, nullptr, -1.0)));
635  }
636  }
637  for (auto const & l : signalNew) {
638  int chan = l.first;
639  auto iter = theSignal.find(chan);
640  if (iter != theSignal.end()) {
641  theSignal[chan] += l.second.ampl();
642  } else {
643  theSignal.insert(std::pair<int,DigitizerUtility::Amplitude>(chan, DigitizerUtility::Amplitude(l.second.ampl(), nullptr, -1.0)));
644  }
645  }
646  }
647  if (!addNoisyPixels) // Option to skip noise in non-hit pixels
648  return;
649 
650  // Add noise on non-hit pixels
651  // Use here the pixel noise
652  int numberOfPixels = numRows * numColumns;
653  std::map<int,float, std::less<int> > otherPixels;
654  std::map<int,float, std::less<int> >::iterator mapI;
655 
656  theNoiser->generate(numberOfPixels,
657  thePixelThreshold, //thr. in un. of nois
658  theNoiseInElectrons, // noise in elec.
659  otherPixels,
660  rengine_ );
661 
662  LogDebug("Phase2TrackerDigitizerAlgorithm")
663  << " Add noisy pixels " << numRows << " "
664  << numColumns << " " << theNoiseInElectrons << " "
665  << theThresholdInE_Endcap << " " << theThresholdInE_Barrel << " " << numberOfPixels << " "
666  << otherPixels.size() ;
667 
668  // Add noisy pixels
669  for (mapI = otherPixels.begin(); mapI!= otherPixels.end(); mapI++) {
670  int iy = ((*mapI).first) / numRows;
671  int ix = ((*mapI).first) - (iy*numRows);
672  // Keep for a while for testing.
673  if( iy < 0 || iy > (numColumns-1) )
674  LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in iy " << iy;
675  if( ix < 0 || ix > (numRows-1) )
676  LogWarning("Phase2TrackerDigitizerAlgorithm") << " error in ix " << ix;
677 
678  int chan;
680 
681  LogDebug ("Phase2TrackerDigitizerAlgorithm")
682  <<" Storing noise = " << (*mapI).first << " " << (*mapI).second
683  << " " << ix << " " << iy << " " << chan ;
684 
685  if (theSignal[chan] == 0) {
686  int noise = int((*mapI).second);
687  theSignal[chan] = DigitizerUtility::Amplitude (noise, nullptr, -1.);
688  }
689  }
690 }
691 
692 // ============================================================================
693 //
694 // Simulate the readout inefficiencies.
695 // Delete a selected number of single pixels, dcols and rocs.
697  const Phase2TrackerGeomDetUnit* pixdet,
698  const TrackerTopology *tTopo) {
699 
700  uint32_t detID = pixdet->geographicalId().rawId();
701 
702  signal_map_type& theSignal = _signal[detID]; // check validity
703 
704  // Predefined efficiencies
705  float subdetEfficiency = 1.0;
706 
707  // setup the chip indices conversion
708  unsigned int Subid=DetId(detID).subdetId();
709  if (Subid == PixelSubdetector::PixelBarrel || Subid == StripSubdetector::TOB) { // barrel layers
710  unsigned int layerIndex = tTopo->pxbLayer(detID);
711  if (layerIndex-1 < eff.barrel_efficiencies.size()) subdetEfficiency = eff.barrel_efficiencies[layerIndex-1];
712  } else { // forward disks
713  unsigned int diskIndex = 2 * tTopo->pxfDisk(detID) - tTopo->pxfSide(detID);
714  if (diskIndex-1 < eff.endcap_efficiencies.size()) subdetEfficiency = eff.endcap_efficiencies[diskIndex-1];
715  }
716 
717  LogDebug ("Phase2TrackerDigitizerAlgorithm") << " enter pixel_inefficiency " << subdetEfficiency;
718 
719  // Now loop again over pixels to kill some of them.
720  // Loop over hits, amplitude in electrons, channel = coded row,col
721  for (auto & s : theSignal) {
722  float rand = flatDistribution_->fire();
723  if( rand>subdetEfficiency ) {
724  // make amplitude =0
725  s.second.set(0.); // reset amplitude,
726  }
727  }
728 }
729 // =======================================================================================
730 //
731 // Set the drift direction accoring to the Bfield in local det-unit frame
732 // Works for both barrel and forward pixels.
733 // Replace the sign convention to fit M.Swartz's formulaes.
734 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
735 // parameter value
736 
738  const GlobalVector& bfield,
739  const DetId& detId) const {
740  Frame detFrame(pixdet->surface().position(),pixdet->surface().rotation());
741  LocalVector Bfield = detFrame.toLocal(bfield);
742  float alpha2_Endcap;
743  float alpha2_Barrel;
744  float alpha2;
745 
746  float dir_x = 0.0;
747  float dir_y = 0.0;
748  float dir_z = 0.0;
749  float scale = 0.0;
750 
751  uint32_t detID = pixdet->geographicalId().rawId();
752  unsigned int Sub_detid = DetId(detID).subdetId();
753 
754  // Read Lorentz angle from cfg file:
755  if (!use_LorentzAngle_DB_) {
756  if (alpha2Order) {
759  } else {
760  alpha2_Endcap = 0.0;
761  alpha2_Barrel = 0.0;
762  }
763 
764  if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) { // barrel layers
765  dir_x = -( tanLorentzAnglePerTesla_Barrel * Bfield.y() + alpha2_Barrel* Bfield.z()* Bfield.x() );
766  dir_y = +( tanLorentzAnglePerTesla_Barrel * Bfield.x() - alpha2_Barrel* Bfield.z()* Bfield.y() );
767  dir_z = -(1 + alpha2_Barrel* Bfield.z()*Bfield.z() );
768  scale = (1 + alpha2_Barrel* Bfield.z()*Bfield.z() );
769 
770  } else { // forward disks
771  dir_x = -( tanLorentzAnglePerTesla_Endcap * Bfield.y() + alpha2_Endcap* Bfield.z()* Bfield.x() );
772  dir_y = +( tanLorentzAnglePerTesla_Endcap * Bfield.x() - alpha2_Endcap* Bfield.z()* Bfield.y() );
773  dir_z = -(1 + alpha2_Endcap* Bfield.z()*Bfield.z() );
774  scale = (1 + alpha2_Endcap* Bfield.z()*Bfield.z() );
775  }
776  }
777 
778  // Read Lorentz angle from DB:
779  if (use_LorentzAngle_DB_) {
780  float lorentzAngle = SiPixelLorentzAngle_->getLorentzAngle(detId);
781  alpha2 = lorentzAngle * lorentzAngle;
782 
783  dir_x = -( lorentzAngle * Bfield.y() + alpha2 * Bfield.z()* Bfield.x() );
784  dir_y = +( lorentzAngle * Bfield.x() - alpha2 * Bfield.z()* Bfield.y() );
785  dir_z = -(1 + alpha2 * Bfield.z()*Bfield.z() );
786  scale = (1 + alpha2 * Bfield.z()*Bfield.z() );
787  }
788 
789  LocalVector theDriftDirection = LocalVector(dir_x/scale, dir_y/scale, dir_z/scale );
790 
791  LogDebug ("Phase2TrackerDigitizerAlgorithm") << " The drift direction in local coordinate is "
792  << theDriftDirection ;
793  return theDriftDirection;
794 }
795 
796 // =============================================================================
797 
799  signal_map_type& theSignal = _signal[detID]; // check validity
800 
801  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
802  for (auto & s : theSignal) {
803  std::pair<int,int> ip;
804  if (pixelFlag) ip = PixelDigi::channelToPixel(s.first);//get pixel pos
805  else ip = Phase2TrackerDigi::channelToPixel(s.first);//get pixel pos
806  int row = ip.first; // X in row
807  int col = ip.second; // Y is in col
808  //transform to ROC index coordinates
809  if (theSiPixelGainCalibrationService_->isDead(detID, col, row)) {
810  s.second.set(0.); // reset amplitude,
811  }
812  }
813 }
814 
815 // ==========================================================================
816 
818  bool isbad = false;
819  int detid = detID;
821  for ( auto const & det_m : DeadModules) {
822  int Dead_detID = det_m.getParameter<int>("Dead_detID");
823  Module = det_m.getParameter<std::string>("Module");
824  if (detid == Dead_detID){
825  isbad = true;
826  break;
827  }
828  }
829 
830  if (!isbad) return;
831 
832  signal_map_type& theSignal = _signal[detID]; // check validity
833 
834 
835  for (auto & s : theSignal) {
836  std::pair<int,int> ip;
837  if (pixelFlag) ip = PixelDigi::channelToPixel(s.first);
838  else ip = Phase2TrackerDigi::channelToPixel(s.first);//get pixel pos
839 
840  if (Module == "whole") s.second.set(0.); // reset amplitude
841  else if (Module == "tbmA" && ip.first >= 80 && ip.first <= 159) s.second.set(0.);
842  else if (Module == "tbmB" && ip.first <= 79) s.second.set(0.);
843  }
844 }
845 // ==========================================================================
847  bool isbad = false;
848 
849  std::vector<SiPixelQuality::disabledModuleType>disabledModules = SiPixelBadModule_->getBadComponentList();
850 
852  for (size_t id = 0;id < disabledModules.size(); id++) {
853  if (detID == disabledModules[id].DetID) {
854  isbad = true;
855  badmodule = disabledModules[id];
856  break;
857  }
858  }
859 
860  if (!isbad)
861  return;
862 
863 
864  signal_map_type& theSignal = _signal[detID]; // check validity
865 
866  if (badmodule.errorType == 0) { // this is a whole dead module.
867  for (auto & s : theSignal) {
868  s.second.set(0.); // reset amplitude
869  }
870  }
871  else { // all other module types: half-modules and single ROCs.
872  // Get Bad ROC position:
873  // follow the example of getBadRocPositions in CondFormats/SiPixelObjects/src/SiPixelQuality.cc
874  std::vector<GlobalPixel> badrocpositions (0);
875  for (unsigned int j = 0; j < 16; j++) {
876  if (SiPixelBadModule_->IsRocBad(detID, j) == true) {
877  std::vector<CablingPathToDetUnit> path = map_.product()->pathToDetUnit(detID);
878  for (auto const & p : path) {
879  const PixelROC* myroc = map_.product()->findItem(p);
880  if (myroc->idInDetUnit() == j) {
881  LocalPixel::RocRowCol local = {39, 25}; //corresponding to center of ROC row, col
882  GlobalPixel global = myroc->toGlobal(LocalPixel(local));
883  badrocpositions.push_back(global);
884  break;
885  }
886  }
887  }
888  }
889 
890 
891  for (auto & s : theSignal) {
892  std::pair<int,int> ip;
893  if (pixelFlag) ip = PixelDigi::channelToPixel(s.first);
894  else ip = Phase2TrackerDigi::channelToPixel(s.first);
895 
896  for (auto const & p : badrocpositions) {
897  for (auto & k : badPixels ) {
898  if ( p.row == k.getParameter<int>("row") &&
899  ip.first == k.getParameter<int>("row") &&
900  std::abs(ip.second - p.col) < k.getParameter<int>("col")) {s.second.set(0.);}
901  }
902  }
903  }
904  }
905 }
907  std::map<int, DigitizerUtility::DigiSimInfo> & digi_map,
908  const TrackerTopology* tTopo)
909 {
910  uint32_t detID = pixdet->geographicalId().rawId();
911  auto it = _signal.find(detID);
912  if (it == _signal.end()) return;
913 
914  const signal_map_type& theSignal = _signal[detID];
915 
916 
917  unsigned int Sub_detid = DetId(detID).subdetId();
918 
919  float theThresholdInE = 0.;
920  float theHIPThresholdInE = 0.;
921  // Define Threshold
922  if (Sub_detid == PixelSubdetector::PixelBarrel || Sub_detid == StripSubdetector::TOB) { // Barrel modules
923  if (addThresholdSmearing) theThresholdInE = smearedThreshold_Barrel_->fire(); // gaussian smearing
924  else theThresholdInE = theThresholdInE_Barrel; // no smearing
925  theHIPThresholdInE = theHIPThresholdInE_Barrel;
926  } else { // Forward disks modules
927  if (addThresholdSmearing) theThresholdInE = smearedThreshold_Endcap_->fire(); // gaussian smearing
928  else theThresholdInE = theThresholdInE_Endcap; // no smearing
929  theHIPThresholdInE = theHIPThresholdInE_Endcap;
930  }
931 
932  if (addNoise) add_noise(pixdet, theThresholdInE/theNoiseInElectrons); // generate noise
933 
934  // Do only if needed
935  if (AddPixelInefficiency && !theSignal.empty()) {
936  if (use_ineff_from_db_)
937  pixel_inefficiency_db(detID);
938  else
939  pixel_inefficiency(subdetEfficiencies_, pixdet, tTopo);
940  }
941  if (use_module_killing_) {
942  if (use_deadmodule_DB_) // remove dead modules using DB
943  module_killing_DB(detID);
944  else // remove dead modules using the list in cfg file
945  module_killing_conf(detID);
946  }
947 
948  // Digitize if the signal is greater than threshold
949  for (auto const & s : theSignal) {
950  // DigitizerUtility::Amplitude sig_data = s.second;
951  const DigitizerUtility::Amplitude& sig_data = s.second;
952  float signalInElectrons = sig_data.ampl();
953  int adc;
954  if (signalInElectrons >= theThresholdInE) { // check threshold
956  else adc = std::min( int(signalInElectrons / theElectronPerADC), theAdcFullScale );
958  info.sig_tot = adc;
959  info.ot_bit = ( signalInElectrons > theHIPThresholdInE ? true : false);
960  if (makeDigiSimLinks_ ) {
961  for (auto const & l : sig_data.simInfoList()) {
962  float charge_frac = l.first/signalInElectrons;
963  if (l.first > -5.0) info.simInfoList.push_back({charge_frac, l.second.get()});
964  }
965  }
966  digi_map.insert({s.first, info});
967  }
968  }
969 }
void primary_ionization(const PSimHit &hit, std::vector< DigitizerUtility::EnergyDepositUnit > &ionization_points) const
int adc(sample_type sample)
get the ADC sample (12 bits)
#define LogDebug(id)
T getParameter(std::string const &) const
virtual int nrows() const =0
LocalVector DriftDirection(const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const DetId &detId) const
static const TGPicture * info(bool iBackgroundIsBlack)
Local3DVector LocalVector
Definition: LocalVector.h:12
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:72
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
T y() const
Definition: PV2DBase.h:46
T perp() const
Definition: PV3DBase.h:72
virtual void add_noise(const Phase2TrackerGeomDetUnit *pixdet, float thePixelThreshold)
unsigned int pxfDisk(const DetId &id) const
void drift(const PSimHit &hit, const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const std::vector< DigitizerUtility::EnergyDepositUnit > &ionization_points, std::vector< DigitizerUtility::SignalPoint > &collection_points) const
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate
T y() const
Definition: PV3DBase.h:63
virtual std::pair< float, float > pitch() const =0
edm::ESHandle< SiPixelFedCablingMap > map_
const Bounds & bounds() const
Definition: Surface.h:120
bool IsRocBad(const uint32_t &detid, const short &rocNb) const
#define nullptr
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
SigmaX
Definition: gun_cff.py:25
virtual void module_killing_conf(uint32_t detID)
identify pixel inside single ROC
Definition: LocalPixel.h:7
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:68
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
uint32_t rawId() const
get the raw id
Definition: DetId.h:44
const std::unique_ptr< SiPixelGainCalibrationOfflineSimService > theSiPixelGainCalibrationService_
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
Phase2TrackerDigitizerAlgorithm(const edm::ParameterSet &conf_common, const edm::ParameterSet &conf_specific, CLHEP::HepRandomEngine &)
T mag() const
Definition: PV3DBase.h:67
std::vector< edm::ParameterSet > badPixels
T sqrt(T t)
Definition: SSEVec.h:18
static PackedDigiType pixelToChannel(unsigned int row, unsigned int col)
edm::ESHandle< SiPixelLorentzAngle > SiPixelLorentzAngle_
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:40
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
virtual int channel(const LocalPoint &p) const =0
virtual void pixel_inefficiency_db(uint32_t detID)
const std::vector< disabledModuleType > getBadComponentList() const
const std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Barrel_
T min(T a, T b)
Definition: MathUtil.h:58
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:38
const 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
ii
Definition: cuy.py:588
#define M_PI
std::vector< edm::ParameterSet > Parameters
int k[5][pyjets_maxn]
unsigned int pxbLayer(const DetId &id) const
void fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments, std::vector< float > &elossVector) const
Definition: DetId.h:18
void induce_signal(const PSimHit &hit, const size_t hitIndex, const unsigned int tofBin, const Phase2TrackerGeomDetUnit *pixdet, const std::vector< DigitizerUtility::SignalPoint > &collection_points)
virtual float thickness() const =0
SigmaY
Definition: gun_cff.py:26
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:62
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
row and collumn in ROC representation
Definition: LocalPixel.h:15
const sipixelobjects::PixelROC * findItem(const sipixelobjects::CablingPathToDetUnit &path) const final
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi...
std::vector< sipixelobjects::CablingPathToDetUnit > pathToDetUnit(uint32_t rawDetId) const final
#define Module(md)
Definition: vmac.h:203
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
HLT enums.
virtual int ncolumns() const =0
int particleType() const
Definition: PSimHit.h:85
col
Definition: cuy.py:1008
Signal rand(Signal arg)
Definition: vlib.cc:442
float getLorentzAngle(const uint32_t &) const
const std::unique_ptr< CLHEP::RandFlat > flatDistribution_
unsigned int pxfSide(const DetId &id) const
const RotationType & rotation() const
const std::unique_ptr< CLHEP::RandGaussQ > gaussDistribution_
std::vector< std::pair< float, SimHitInfoForLinks * > > simInfoList
virtual void pixel_inefficiency(const SubdetEfficiencies &eff, const Phase2TrackerGeomDetUnit *pixdet, const TrackerTopology *tTopo)
const std::vector< std::pair< float, std::unique_ptr< SimHitInfoForLinks > > > & simInfoList() const
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser
T x() const
Definition: PV2DBase.h:45
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
T const * product() const
Definition: ESHandle.h:86
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
static std::pair< unsigned int, unsigned int > channelToPixel(unsigned int ch)
*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:93
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:45
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:59
edm::ESHandle< SiPixelQuality > SiPixelBadModule_