CMS 3D CMS Logo

SiPixelDigitizerAlgorithm.h
Go to the documentation of this file.
1 #ifndef SiPixelDigitizerAlgorithm_h
2 #define SiPixelDigitizerAlgorithm_h
3 
4 #include <map>
5 #include <memory>
6 #include <vector>
7 #include <iostream>
30 #include "boost/multi_array.hpp"
31 
32 typedef boost::multi_array<float, 2> array_2d;
33 
34 // forward declarations
35 
36 // For the random numbers
37 namespace CLHEP {
38  class HepRandomEngine;
39 }
40 
41 class DetId;
43 class PixelDigi;
44 class PixelDigiSimLink;
45 class PixelGeomDetUnit;
50 class SiPixelQuality;
52 class TrackerGeometry;
53 class TrackerTopology;
57 
59 public:
62 
63  // initialization that cannot be done in the constructor
64  void init(const edm::EventSetup& es);
65 
66  void initializeEvent() { _signal.clear(); }
67 
68  //run the algorithm to digitize a single det
69  void accumulateSimHits(const std::vector<PSimHit>::const_iterator inputBegin,
70  const std::vector<PSimHit>::const_iterator inputEnd,
71  const size_t inputBeginGlobalIndex,
72  const unsigned int tofBin,
73  const PixelGeomDetUnit* pixdet,
74  const GlobalVector& bfield,
75  const TrackerTopology* tTopo,
76  CLHEP::HepRandomEngine*);
77  void digitize(const PixelGeomDetUnit* pixdet,
78  std::vector<PixelDigi>& digis,
79  std::vector<PixelDigiSimLink>& simlinks,
80  std::vector<PixelDigiAddTempInfo>& newClass_Digi_extra,
81  const TrackerTopology* tTopo,
82  CLHEP::HepRandomEngine*);
84  void fillSimHitMaps(std::vector<PSimHit> simHits, const unsigned int tofBin);
85  void resetSimHitMaps();
86  void init_DynIneffDB(const edm::EventSetup&);
87  std::unique_ptr<PixelFEDChannelCollection> chooseScenario(PileupMixingContent* puInfo, CLHEP::HepRandomEngine*);
88 
89  void lateSignalReweight(const PixelGeomDetUnit* pixdet,
90  std::vector<PixelDigi>& digis,
91  std::vector<PixelSimHitExtraInfo>& newClass_Sim_extra,
92  const TrackerTopology* tTopo,
93  CLHEP::HepRandomEngine* engine);
94 
95  // for premixing
96  void calculateInstlumiFactor(const std::vector<PileupSummaryInfo>& ps,
97  int bunchSpacing); // TODO: try to remove the duplication of logic...
98  void setSimAccumulator(const std::map<uint32_t, std::map<int, int> >& signalMap);
99  std::unique_ptr<PixelFEDChannelCollection> chooseScenario(const std::vector<PileupSummaryInfo>& ps,
100  CLHEP::HepRandomEngine* engine);
101 
102  bool killBadFEDChannels() const;
103  typedef std::unordered_map<std::string, PixelFEDChannelCollection> PixelFEDChannelCollectionMap;
105 
106  class Amplitude {
107  public:
108  Amplitude() : _amp(0.0) {}
109  Amplitude(float amp, float frac) : _amp(amp), _frac(1, frac) {
110  //in case of digi from noisypixels
111  //the MC information are removed
112  if (_frac[0] < -0.5) {
113  _frac.pop_back();
114  }
115  }
116 
117  Amplitude(float amp, const PSimHit* hitp, size_t hitIndex, size_t hitInd4CR, unsigned int tofBin, float frac)
118  : _amp(amp), _frac(1, frac) {
119  //in case of digi from noisypixels
120  //the MC information are removed
121  if (_frac[0] < -0.5) {
122  _frac.pop_back();
123  } else {
124  _hitInfos.emplace_back(hitp, hitIndex, tofBin, hitInd4CR, amp);
125  }
126  }
127 
128  // can be used as a float by convers.
129  operator float() const { return _amp; }
130  float ampl() const { return _amp; }
131  const std::vector<float>& individualampl() const { return _frac; }
132  const std::vector<SimHitInfoForLinks>& hitInfos() const { return _hitInfos; }
133 
134  void operator+=(const Amplitude& other) {
135  _amp += other._amp;
136  //in case of contribution of noise to the digi
137  //the MC information are removed
138  if (other._frac[0] > -0.5) {
139  if (!other._hitInfos.empty()) {
140  _hitInfos.insert(_hitInfos.end(), other._hitInfos.begin(), other._hitInfos.end());
141  }
142  _frac.insert(_frac.end(), other._frac.begin(), other._frac.end());
143  }
144  }
145  void operator+=(const float& amp) { _amp += amp; }
146 
147  void set(const float amplitude) { // Used to reset the amplitude
148  _amp = amplitude;
149  }
150  /* void setind (const float indamplitude) { // Used to reset the amplitude */
151  /* _frac = idamplitude; */
152  /* } */
153  private:
154  float _amp;
155  std::vector<float> _frac;
156  std::vector<SimHitInfoForLinks> _hitInfos;
157  }; // end class Amplitude
158 
159 private:
160  //Accessing Lorentz angle from DB:
163 
164  //Accessing Dead pixel modules from DB:
167 
168  //Accessing Map and Geom:
171  const SiPixelFedCablingMap* map_ = nullptr;
172  const TrackerGeometry* geom_ = nullptr;
173 
174  // Get Dynamic Inefficiency scale factors from DB
177 
178  // For BadFEDChannel simulation
183  // Define internal classes
184 
185  // definition class
186  //
187 
188  // Define a class to hold the calibration parameters per pixel
189  // Internal
191  public:
192  float p0;
193  float p1;
194  float p2;
195  float p3;
196  };
197  //
198  // Define a class for 3D ionization points and energy
199  //
204  public:
205  EnergyDepositUnit() : _energy(0), _position(0, 0, 0) {}
206  EnergyDepositUnit(float energy, float x, float y, float z) : _energy(energy), _position(x, y, z) {}
208  float x() const { return _position.x(); }
209  float y() const { return _position.y(); }
210  float z() const { return _position.z(); }
211  float energy() const { return _energy; }
212 
213  private:
214  float _energy;
216  };
217 
218  //
219  // define class to store signals on the collection surface
220  //
225  class SignalPoint {
226  public:
227  SignalPoint() : _pos(0, 0), _time(0), _amplitude(0), _sigma_x(1.), _sigma_y(1.), _hitp(nullptr) {}
228 
229  SignalPoint(float x, float y, float sigma_x, float sigma_y, float t, float a = 1.0)
230  : _pos(x, y), _time(t), _amplitude(a), _sigma_x(sigma_x), _sigma_y(sigma_y), _hitp(nullptr) {}
231 
232  SignalPoint(float x, float y, float sigma_x, float sigma_y, float t, const PSimHit& hit, float a = 1.0)
234 
235  const LocalPoint& position() const { return _pos; }
236  float x() const { return _pos.x(); }
237  float y() const { return _pos.y(); }
238  float sigma_x() const { return _sigma_x; }
239  float sigma_y() const { return _sigma_y; }
240  float time() const { return _time; }
241  float amplitude() const { return _amplitude; }
242  const PSimHit& hit() { return *_hitp; }
244  _amplitude = amp;
245  return *this;
246  }
247 
248  private:
250  float _time;
251  float _amplitude;
252  float _sigma_x; // gaussian sigma in the x direction (cm)
253  float _sigma_y; // " " y direction (cm) */
254  const PSimHit* _hitp;
255  };
256 
257  //
258  // PixelEfficiencies struct
259  //
267  int NumberOfEndcapDisks);
268  bool FromConfig; // If true read from Config, otherwise use Database
269 
271  std::vector<double> pu_scale; // in config: 0-3 BPix, 4-5 FPix (inner, outer)
272  std::vector<std::vector<double> > thePUEfficiency; // Instlumi dependent efficiency
273 
274  // Read factors from Configuration
275  double thePixelEfficiency[20]; // Single pixel effciency
276  double thePixelColEfficiency[20]; // Column effciency
277  double thePixelChipEfficiency[20]; // ROC efficiency
278  std::vector<double> theLadderEfficiency_BPix[20]; // Ladder efficiency
279  std::vector<double> theModuleEfficiency_BPix[20]; // Module efficiency
280  double theInnerEfficiency_FPix[20]; // Fpix inner module efficiency
281  double theOuterEfficiency_FPix[20]; // Fpix outer module efficiency
282  unsigned int FPixIndex; // The Efficiency index for FPix Disks
283 
284  // Read factors from DB and fill containers
285  std::map<uint32_t, double> PixelGeomFactors;
286  std::map<uint32_t, std::vector<double> > PixelGeomFactorsROCStdPixels;
287  std::map<uint32_t, std::vector<double> > PixelGeomFactorsROCBigPixels;
288  std::map<uint32_t, double> ColGeomFactors;
289  std::map<uint32_t, double> ChipGeomFactors;
290  std::map<uint32_t, size_t> iPU;
291 
292  // constants for ROC level simulation for Phase1
294  static const int rocIdMaskBits = 0x1F;
296  bool matches(const DetId&, const DetId&, const std::vector<uint32_t>&);
297  std::unique_ptr<PixelFEDChannelCollection> PixelFEDChannelCollection_;
298  };
299 
300  //
301  // PixelAging struct
302  //
306  struct PixelAging {
308  float thePixelPseudoRadDamage[20]; // PseudoRadiation Damage Values for aging studies
309  unsigned int FPixIndex; // The Efficiency index for FPix Disks
310  };
311 
312 private:
313  // Internal typedefs
314  typedef std::map<int, Amplitude, std::less<int> > signal_map_type; // from Digi.Skel.
315  typedef signal_map_type::iterator signal_map_iterator; // from Digi.Skel.
316  typedef signal_map_type::const_iterator signal_map_const_iterator; // from Digi.Skel.
317  typedef std::map<uint32_t, signal_map_type> signalMaps;
319  typedef std::vector<edm::ParameterSet> Parameters;
320  typedef boost::multi_array<float, 2> array_2d;
321 
322  typedef std::pair<unsigned int, unsigned int> subDetTofBin;
323  typedef std::map<unsigned int, std::vector<PSimHit> > simhit_map;
325  typedef std::map<subDetTofBin, unsigned int> simhit_collectionMap;
327 
328  // Contains the accumulated hit info.
330 
331  const bool makeDigiSimLinks_;
333 
334  const bool use_ineff_from_db_;
335  const bool use_module_killing_; // remove or not the dead pixel modules
336  const bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
337  const bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
338 
340 
341  std::unique_ptr<SiPixelChargeReweightingAlgorithm> TheNewSiPixelChargeReweightingAlgorithmClass;
342 
343 private:
344  // Variables
345  //external parameters
346  // go from Geant energy GeV to number of electrons
347  const float GeVperElectron; // 3.7E-09
348 
349  //-- drift
350  const float Sigma0; //=0.0007 // Charge diffusion in microns for 300 micron Si
351  const float Dist300; //=0.0300 // Define 300microns for normalization
352  const bool alpha2Order; // Switch on/off of E.B effect
353 
354  //-- induce_signal
355  const float ClusterWidth; // Gaussian charge cutoff width in sigma units
356  //-- Allow for upgrades
357  const int NumberOfBarrelLayers; // Default = 3
358  const int NumberOfEndcapDisks; // Default = 2
359 
360  //-- make_digis
361  const float theElectronPerADC; // Gain, number of electrons per adc count.
362  const int theAdcFullScale; // Saturation count, 255=8bit.
363  const int theAdcFullScLateCR; // Saturation count, 255=8bit.
364  const float theNoiseInElectrons; // Noise (RMS) in units of electrons.
365  const float theReadoutNoise; // Noise of the readount chain in elec,
366  //inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
367 
368  const float theThresholdInE_FPix; // Pixel threshold in electrons FPix.
369  const float theThresholdInE_BPix; // Pixel threshold in electrons BPix.
370  const float theThresholdInE_BPix_L1; // In case the BPix layer1 gets a different threshold
371  const float theThresholdInE_BPix_L2; // In case the BPix layer2 gets a different threshold
372 
377 
378  const float electronsPerVCAL; // for electrons - VCAL conversion
379  const float electronsPerVCAL_Offset; // in misscalibrate()
380  const float electronsPerVCAL_L1; // same for Layer 1
381  const float electronsPerVCAL_L1_Offset; // same for Layer 1
382 
383  const float theTofLowerCut; // Cut on the particle TOF
384  const float theTofUpperCut; // Cut on the particle TOF
385  const float tanLorentzAnglePerTesla_FPix; //FPix Lorentz angle tangent per Tesla
386  const float tanLorentzAnglePerTesla_BPix; //BPix Lorentz angle tangent per Tesla
387 
388  const float FPix_p0;
389  const float FPix_p1;
390  const float FPix_p2;
391  const float FPix_p3;
392  const float BPix_p0;
393  const float BPix_p1;
394  const float BPix_p2;
395  const float BPix_p3;
396 
397  //-- add_noise
398  const bool addNoise;
400  const bool addNoisyPixels;
401  const bool fluctuateCharge;
402 
403  //-- pixel efficiency
404  const bool AddPixelInefficiency; // bool to read in inefficiencies
405  const bool KillBadFEDChannels;
407 
408  //-- calibration smearing
409  const bool doMissCalibrate; // Switch on the calibration smearing
410  const bool doMissCalInLateCR; // Switch on the calibration smearing
411  const float theGainSmearing; // The sigma of the gain fluctuation (around 1)
412  const float theOffsetSmearing; // The sigma of the offset fluct. (around 0)
413 
414  // pixel aging
415  const bool AddPixelAging;
416  const bool UseReweighting;
417 
418  // The PDTable
419  //HepPDTable *particleTable;
420  //ParticleDataTable *particleTable;
421 
422  //-- charge fluctuation
423  const double tMax; // The delta production cut, should be as in OSCAR = 30keV
424  // cmsim = 100keV
425 
426  // The eloss fluctuation class from G4. Is the right place?
427  const std::unique_ptr<SiG4UniversalFluctuation> fluctuate; // make a pointer
428  const std::unique_ptr<GaussianTailNoiseGenerator> theNoiser;
429 
430  // To store calibration constants
431  const std::map<int, CalParameters, std::less<int> > calmap;
432 
433  //-- additional member functions
434  // Private methods
435  std::map<int, CalParameters, std::less<int> > initCal() const;
436  void primary_ionization(const PSimHit& hit,
437  std::vector<EnergyDepositUnit>& ionization_points,
438  CLHEP::HepRandomEngine*) const;
439  void drift(const PSimHit& hit,
440  const PixelGeomDetUnit* pixdet,
441  const GlobalVector& bfield,
442  const TrackerTopology* tTopo,
443  const std::vector<EnergyDepositUnit>& ionization_points,
444  std::vector<SignalPoint>& collection_points) const;
445  void induce_signal(std::vector<PSimHit>::const_iterator inputBegin,
446  std::vector<PSimHit>::const_iterator inputEnd,
447  const PSimHit& hit,
448  const size_t hitIndex,
449  const size_t FirstHitIndex,
450  const unsigned int tofBin,
451  const PixelGeomDetUnit* pixdet,
452  const std::vector<SignalPoint>& collection_points);
453  void fluctuateEloss(int particleId,
454  float momentum,
455  float eloss,
456  float length,
457  int NumberOfSegments,
458  float elossVector[],
459  CLHEP::HepRandomEngine*) const;
460  void add_noise(const PixelGeomDetUnit* pixdet, float thePixelThreshold, CLHEP::HepRandomEngine*);
461  void make_digis(float thePixelThresholdInE,
462  uint32_t detID,
463  const PixelGeomDetUnit* pixdet,
464  std::vector<PixelDigi>& digis,
465  std::vector<PixelDigiSimLink>& simlinks,
466  std::vector<PixelDigiAddTempInfo>& newClass_Digi_extra,
467  const TrackerTopology* tTopo) const;
468  void pixel_inefficiency(const PixelEfficiencies& eff,
469  const PixelGeomDetUnit* pixdet,
470  const TrackerTopology* tTopo,
471  CLHEP::HepRandomEngine*);
472 
473  void pixel_inefficiency_db(uint32_t detID);
474 
475  float pixel_aging(const PixelAging& aging, const PixelGeomDetUnit* pixdet, const TrackerTopology* tTopo) const;
476 
477  // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
478  const std::unique_ptr<SiPixelGainCalibrationOfflineSimService> theSiPixelGainCalibrationService_;
479  float missCalibrate(
480  uint32_t detID, const TrackerTopology* tTopo, const PixelGeomDetUnit* pixdet, int col, int row, float amp) const;
481  LocalVector DriftDirection(const PixelGeomDetUnit* pixdet, const GlobalVector& bfield, const DetId& detId) const;
482 
483  void module_killing_conf(
484  uint32_t detID); // remove dead modules using the list in the configuration file PixelDigi_cfi.py
485  void module_killing_DB(uint32_t detID); // remove dead modules uisng the list in the DB
486 
489 
490  double calcQ(float x) const {
491  // need erf(x/sqrt2)
492  //float x2=0.5*x*x;
493  //float a=0.147;
494  //double erf=sqrt(1.0f-exp( -1.0f*x2*( (4/M_PI)+a*x2)/(1.0+a*x2)));
495  //if (x<0.) erf*=-1.0;
496  //return 0.5*(1.0-erf);
497 
498  auto xx = std::min(0.5f * x * x, 12.5f);
499  return 0.5 * (1.0 - std::copysign(std::sqrt(1.f - unsafe_expf<4>(-xx * (1.f + 0.2733f / (1.f + 0.147f * xx)))), x));
500  }
501 };
502 
503 #endif
void init(const edm::EventSetup &es)
void fillSimHitMaps(std::vector< PSimHit > simHits, const unsigned int tofBin)
void pixel_inefficiency_db(uint32_t detID)
signal_map_type::const_iterator signal_map_const_iterator
SignalPoint(float x, float y, float sigma_x, float sigma_y, float t, const PSimHit &hit, float a=1.0)
const std::unique_ptr< SiPixelGainCalibrationOfflineSimService > theSiPixelGainCalibrationService_
EnergyDepositUnit(float energy, Local3DPoint position)
T z() const
Definition: PV3DBase.h:61
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate
PixelEfficiencies(const edm::ParameterSet &conf, bool AddPixelInefficiency, int NumberOfBarrelLayers, int NumberOfEndcapDisks)
SignalPoint(float x, float y, float sigma_x, float sigma_y, float t, float a=1.0)
edm::ESGetToken< SiPixelQualityProbabilities, SiPixelStatusScenarioProbabilityRcd > scenarioProbabilityToken_
std::unique_ptr< PixelFEDChannelCollection > PixelFEDChannelCollection_
std::vector< std::vector< double > > thePUEfficiency
void make_digis(float thePixelThresholdInE, uint32_t detID, const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, std::vector< PixelDigiSimLink > &simlinks, std::vector< PixelDigiAddTempInfo > &newClass_Digi_extra, const TrackerTopology *tTopo) const
float pixel_aging(const PixelAging &aging, const PixelGeomDetUnit *pixdet, const TrackerTopology *tTopo) const
void drift(const PSimHit &hit, const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const TrackerTopology *tTopo, const std::vector< EnergyDepositUnit > &ionization_points, std::vector< SignalPoint > &collection_points) const
std::map< unsigned int, std::vector< PSimHit > > simhit_map
const SiPixelDynamicInefficiency * SiPixelDynamicInefficiency_
std::unique_ptr< SiPixelChargeReweightingAlgorithm > TheNewSiPixelChargeReweightingAlgorithmClass
std::map< subDetTofBin, unsigned int > simhit_collectionMap
void induce_signal(std::vector< PSimHit >::const_iterator inputBegin, std::vector< PSimHit >::const_iterator inputEnd, const PSimHit &hit, const size_t hitIndex, const size_t FirstHitIndex, const unsigned int tofBin, const PixelGeomDetUnit *pixdet, const std::vector< SignalPoint > &collection_points)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
LocalVector DriftDirection(const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const DetId &detId) const
edm::ESGetToken< SiPixelQuality, SiPixelQualityRcd > SiPixelBadModuleToken_
SiPixelDigitizerAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
const PixelFEDChannelCollectionMap * quality_map
std::pair< unsigned int, unsigned int > subDetTofBin
void digitize(const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, std::vector< PixelDigiSimLink > &simlinks, std::vector< PixelDigiAddTempInfo > &newClass_Digi_extra, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
boost::multi_array< float, 2 > array_2d
std::map< uint32_t, std::vector< double > > PixelGeomFactorsROCBigPixels
void primary_ionization(const PSimHit &hit, std::vector< EnergyDepositUnit > &ionization_points, CLHEP::HepRandomEngine *) const
T sqrt(T t)
Definition: SSEVec.h:19
const SiPixelLorentzAngle * SiPixelLorentzAngle_
double f[11][100]
float missCalibrate(uint32_t detID, const TrackerTopology *tTopo, const PixelGeomDetUnit *pixdet, int col, int row, float amp) const
void init_from_db(const TrackerGeometry *, const SiPixelDynamicInefficiency *)
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleSimRcd > SiPixelLorentzAngleToken_
signal_map_type::iterator signal_map_iterator
void init_DynIneffDB(const edm::EventSetup &)
std::unique_ptr< PixelFEDChannelCollection > chooseScenario(PileupMixingContent *puInfo, CLHEP::HepRandomEngine *)
EnergyDepositUnit(float energy, float x, float y, float z)
boost::multi_array< float, 2 > array_2d
void setSimAccumulator(const std::map< uint32_t, std::map< int, int > > &signalMap)
Amplitude(float amp, const PSimHit *hitp, size_t hitIndex, size_t hitInd4CR, unsigned int tofBin, float frac)
edm::ESGetToken< PixelFEDChannelCollectionMap, SiPixelFEDChannelContainerESProducerRcd > PixelFEDChannelCollectionMapToken_
const std::map< int, CalParameters, std::less< int > > calmap
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser
Definition: DetId.h:17
std::map< int, Amplitude, std::less< int > > signal_map_type
const SiPixelQuality * SiPixelBadModule_
const SiPixelFedCablingMap * map_
std::map< uint32_t, std::vector< double > > PixelGeomFactorsROCStdPixels
void calculateInstlumiFactor(PileupMixingContent *puInfo)
std::unordered_map< std::string, PixelFEDChannelCollection > PixelFEDChannelCollectionMap
void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, const size_t inputBeginGlobalIndex, const unsigned int tofBin, const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
GloballyPositioned< double > Frame
const SiPixelQualityProbabilities * scenarioProbability_
void fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments, float elossVector[], CLHEP::HepRandomEngine *) const
bool matches(const DetId &, const DetId &, const std::vector< uint32_t > &)
void lateSignalReweight(const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, std::vector< PixelSimHitExtraInfo > &newClass_Sim_extra, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *engine)
std::vector< edm::ParameterSet > Parameters
const std::vector< SimHitInfoForLinks > & hitInfos() const
std::vector< SimHitInfoForLinks > _hitInfos
double a
Definition: hdecay.h:119
static int position[264][3]
Definition: ReadPGInfo.cc:289
col
Definition: cuy.py:1009
edm::ESGetToken< SiPixelDynamicInefficiency, SiPixelDynamicInefficiencyRcd > SiPixelDynamicInefficiencyToken_
std::map< int, CalParameters, std::less< int > > initCal() const
PixelAging(const edm::ParameterSet &conf, bool AddPixelAging, int NumberOfBarrelLayers, int NumberOfEndcapDisks)
const edm::ESGetToken< SiPixelFedCablingMap, SiPixelFedCablingMapRcd > mapToken_
std::map< uint32_t, signal_map_type > signalMaps
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
Definition: aging.py:1
const std::vector< float > & individualampl() const
void add_noise(const PixelGeomDetUnit *pixdet, float thePixelThreshold, CLHEP::HepRandomEngine *)
void pixel_inefficiency(const PixelEfficiencies &eff, const PixelGeomDetUnit *pixdet, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)