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