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>
16 
17 // forward declarations
18 
19 
20 // For the random numbers
21 namespace CLHEP {
22  class HepRandomEngine;
23 }
24 
25 namespace edm {
26  class EventSetup;
27  class ParameterSet;
28 }
29 
30 class DetId;
32 class PixelDigi;
33 class PixelDigiSimLink;
34 class PixelGeomDetUnit;
39 class SiPixelQuality;
41 class TrackerGeometry;
42 class TrackerTopology;
43 
45  public:
48 
49  // initialization that cannot be done in the constructor
50  void init(const edm::EventSetup& es);
51 
52  void initializeEvent() {
53  _signal.clear();
54  }
55 
56  //run the algorithm to digitize a single det
57  void accumulateSimHits(const std::vector<PSimHit>::const_iterator inputBegin,
58  const std::vector<PSimHit>::const_iterator inputEnd,
59  const size_t inputBeginGlobalIndex,
60  const unsigned int tofBin,
61  const PixelGeomDetUnit *pixdet,
62  const GlobalVector& bfield,
63  const TrackerTopology *tTopo,
64  CLHEP::HepRandomEngine*);
65  void digitize(const PixelGeomDetUnit *pixdet,
66  std::vector<PixelDigi>& digis,
67  std::vector<PixelDigiSimLink>& simlinks,
68  const TrackerTopology *tTopo,
69  CLHEP::HepRandomEngine*);
70  void calculateInstlumiFactor(PileupMixingContent* puInfo);
71  void init_DynIneffDB(const edm::EventSetup&, const unsigned int&);
72 
73  private:
74 
75  //Accessing Lorentz angle from DB:
77 
78  //Accessing Dead pixel modules from DB:
80 
81  //Accessing Map and Geom:
84 
85  // Get Dynamic Inefficiency scale factors from DB
87 
88  // Define internal classes
89 
90  // definition class
91  //
92  class Amplitude {
93  public:
94  Amplitude() : _amp(0.0) {}
95  Amplitude( float amp, float frac) :
96  _amp(amp), _frac(1, frac) {
97  //in case of digi from noisypixels
98  //the MC information are removed
99  if (_frac[0]<-0.5) {
100  _frac.pop_back();
101  }
102  }
103 
104  Amplitude( float amp, const PSimHit* hitp, size_t hitIndex, unsigned int tofBin, float frac) :
105  _amp(amp), _frac(1, frac) {
106 
107  //in case of digi from noisypixels
108  //the MC information are removed
109  if (_frac[0]<-0.5) {
110  _frac.pop_back();
111  }
112  else {
113  _hitInfos.emplace_back(hitp, hitIndex, tofBin);
114  }
115  }
116 
117  // can be used as a float by convers.
118  operator float() const { return _amp;}
119  float ampl() const {return _amp;}
120  const std::vector<float>& individualampl() const {return _frac;}
121  const std::vector<SimHitInfoForLinks>& hitInfos() const { return _hitInfos; }
122 
123  void operator+=( const Amplitude& other) {
124  _amp += other._amp;
125  //in case of contribution of noise to the digi
126  //the MC information are removed
127  if (other._frac[0]>-0.5){
128  if(!other._hitInfos.empty()) {
129  _hitInfos.insert(_hitInfos.end(), other._hitInfos.begin(), other._hitInfos.end());
130  }
131  _frac.insert(_frac.end(), other._frac.begin(), other._frac.end());
132  }
133  }
134  void operator+=( const float& amp) {
135  _amp += amp;
136  }
137 
138  void set (const float amplitude) { // Used to reset the amplitude
139  _amp = amplitude;
140  }
141 /* void setind (const float indamplitude) { // Used to reset the amplitude */
142 /* _frac = idamplitude; */
143 /* } */
144  private:
145  float _amp;
146  std::vector<float> _frac;
147  std::vector<SimHitInfoForLinks> _hitInfos;
148  }; // end class Amplitude
149 
150  // Define a class to hold the calibration parameters per pixel
151  // Internal
153  public:
154  float p0;
155  float p1;
156  float p2;
157  float p3;
158  };
159  //
160  // Define a class for 3D ionization points and energy
161  //
166  public:
167  EnergyDepositUnit(): _energy(0),_position(0,0,0){}
168  EnergyDepositUnit(float energy,float x, float y, float z):
169  _energy(energy),_position(x,y,z){}
171  _energy(energy),_position(position){}
172  float x() const{return _position.x();}
173  float y() const{return _position.y();}
174  float z() const{return _position.z();}
175  float energy() const { return _energy;}
176  private:
177  float _energy;
179  };
180 
181  //
182  // define class to store signals on the collection surface
183  //
188  class SignalPoint {
189  public:
190  SignalPoint() : _pos(0,0), _time(0), _amplitude(0),
191  _sigma_x(1.), _sigma_y(1.), _hitp(0) {}
192 
193  SignalPoint( float x, float y, float sigma_x, float sigma_y,
194  float t, float a=1.0) :
195  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
196  _sigma_y(sigma_y), _hitp(0) {}
197 
198  SignalPoint( float x, float y, float sigma_x, float sigma_y,
199  float t, const PSimHit& hit, float a=1.0) :
200  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
201  _sigma_y(sigma_y),_hitp(&hit) {}
202 
203  const LocalPoint& position() const { return _pos;}
204  float x() const { return _pos.x();}
205  float y() const { return _pos.y();}
206  float sigma_x() const { return _sigma_x;}
207  float sigma_y() const { return _sigma_y;}
208  float time() const { return _time;}
209  float amplitude() const { return _amplitude;}
210  const PSimHit& hit() { return *_hitp;}
211  SignalPoint& set_amplitude( float amp) { _amplitude = amp; return *this;}
212 
213 
214 
215  private:
217  float _time;
218  float _amplitude;
219  float _sigma_x; // gaussian sigma in the x direction (cm)
220  float _sigma_y; // " " y direction (cm) */
221  const PSimHit* _hitp;
222  };
223 
224  //
225  // PixelEfficiencies struct
226  //
231  PixelEfficiencies(const edm::ParameterSet& conf, bool AddPixelInefficiency, int NumberOfBarrelLayers, int NumberOfEndcapDisks);
232  bool FromConfig; // If true read from Config, otherwise use Database
233 
235  std::vector<double> pu_scale; // in config: 0-3 BPix, 4-5 FPix (inner, outer)
236  std::vector<std::vector<double> > thePUEfficiency; // Instlumi dependent efficiency
237 
238  // Read factors from Configuration
239  double thePixelEfficiency[20]; // Single pixel effciency
240  double thePixelColEfficiency[20]; // Column effciency
241  double thePixelChipEfficiency[20]; // ROC efficiency
242  std::vector<double> theLadderEfficiency_BPix[20]; // Ladder efficiency
243  std::vector<double> theModuleEfficiency_BPix[20]; // Module efficiency
244  double theInnerEfficiency_FPix[20]; // Fpix inner module efficiency
245  double theOuterEfficiency_FPix[20]; // Fpix outer module efficiency
246  unsigned int FPixIndex; // The Efficiency index for FPix Disks
247 
248  // Read factors from DB and fill containers
249  std::map<uint32_t, double> PixelGeomFactors;
250  std::map<uint32_t, double> ColGeomFactors;
251  std::map<uint32_t, double> ChipGeomFactors;
252  std::map<uint32_t, size_t > iPU;
253 
255  bool matches(const DetId&, const DetId&, const std::vector<uint32_t >&);
256  };
257 
258  //
259  // PixelAging struct
260  //
264  struct PixelAging {
265  PixelAging(const edm::ParameterSet& conf, bool AddPixelAging, int NumberOfBarrelLayers, int NumberOfEndcapDisks);
266  float thePixelPseudoRadDamage[20]; // PseudoRadiation Damage Values for aging studies
267  unsigned int FPixIndex; // The Efficiency index for FPix Disks
268  };
269 
270  private:
271  // Internal typedefs
272  typedef std::map<int, Amplitude, std::less<int> > signal_map_type; // from Digi.Skel.
273  typedef signal_map_type::iterator signal_map_iterator; // from Digi.Skel.
274  typedef signal_map_type::const_iterator signal_map_const_iterator; // from Digi.Skel.
275  typedef std::map<uint32_t, signal_map_type> signalMaps;
277  typedef std::vector<edm::ParameterSet> Parameters;
278 
279  // Contains the accumulated hit info.
280  signalMaps _signal;
281 
282  const bool makeDigiSimLinks_;
283 
284  const bool use_ineff_from_db_;
285  const bool use_module_killing_; // remove or not the dead pixel modules
286  const bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
287  const bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
288 
289  const Parameters DeadModules;
290 
291  // Variables
292  //external parameters
293  // go from Geant energy GeV to number of electrons
294  const float GeVperElectron; // 3.7E-09
295 
296  //-- drift
297  const float Sigma0; //=0.0007 // Charge diffusion in microns for 300 micron Si
298  const float Dist300; //=0.0300 // Define 300microns for normalization
299  const bool alpha2Order; // Switch on/off of E.B effect
300 
301 
302 
303  //-- induce_signal
304  const float ClusterWidth; // Gaussian charge cutoff width in sigma units
305  //-- Allow for upgrades
306  const int NumberOfBarrelLayers; // Default = 3
307  const int NumberOfEndcapDisks; // Default = 2
308 
309  //-- make_digis
310  const float theElectronPerADC; // Gain, number of electrons per adc count.
311  const int theAdcFullScale; // Saturation count, 255=8bit.
312  const int theAdcFullScaleStack; // Saturation count for stack layers, 1=1bit.
313  const float theNoiseInElectrons; // Noise (RMS) in units of electrons.
314  const float theReadoutNoise; // Noise of the readount chain in elec,
315  //inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
316 
317  const float theThresholdInE_FPix; // Pixel threshold in electrons FPix.
318  const float theThresholdInE_BPix; // Pixel threshold in electrons BPix.
319  const float theThresholdInE_BPix_L1; // In case the BPix layer1 gets a different threshold
320 
324 
325  const double electronsPerVCAL; // for electrons - VCAL conversion
326  const double electronsPerVCAL_Offset; // in misscalibrate()
327 
328  const float theTofLowerCut; // Cut on the particle TOF
329  const float theTofUpperCut; // Cut on the particle TOF
330  const float tanLorentzAnglePerTesla_FPix; //FPix Lorentz angle tangent per Tesla
331  const float tanLorentzAnglePerTesla_BPix; //BPix Lorentz angle tangent per Tesla
332 
333  const float FPix_p0;
334  const float FPix_p1;
335  const float FPix_p2;
336  const float FPix_p3;
337  const float BPix_p0;
338  const float BPix_p1;
339  const float BPix_p2;
340  const float BPix_p3;
341 
342 
343  //-- add_noise
344  const bool addNoise;
346  const bool addNoisyPixels;
347  const bool fluctuateCharge;
348  //-- pixel efficiency
349  const bool AddPixelInefficiency; // bool to read in inefficiencies
350 
352 
353  //-- calibration smearing
354  const bool doMissCalibrate; // Switch on the calibration smearing
355  const float theGainSmearing; // The sigma of the gain fluctuation (around 1)
356  const float theOffsetSmearing; // The sigma of the offset fluct. (around 0)
357 
358  // pixel aging
359  const bool AddPixelAging;
360 
361  // The PDTable
362  //HepPDTable *particleTable;
363  //ParticleDataTable *particleTable;
364 
365  //-- charge fluctuation
366  const double tMax; // The delta production cut, should be as in OSCAR = 30keV
367  // cmsim = 100keV
368 
369  // The eloss fluctuation class from G4. Is the right place?
370  const std::unique_ptr<SiG4UniversalFluctuation> fluctuate; // make a pointer
371  const std::unique_ptr<GaussianTailNoiseGenerator> theNoiser;
372 
373  // To store calibration constants
374  const std::map<int,CalParameters,std::less<int> > calmap;
375 
376 
377  //-- additional member functions
378  // Private methods
379  std::map<int,CalParameters,std::less<int> > initCal() const;
380  void primary_ionization( const PSimHit& hit, std::vector<EnergyDepositUnit>& ionization_points, CLHEP::HepRandomEngine*) const;
381  void drift(const PSimHit& hit,
382  const PixelGeomDetUnit *pixdet,
383  const GlobalVector& bfield,
384  const TrackerTopology *tTopo,
385  const std::vector<EnergyDepositUnit>& ionization_points,
386  std::vector<SignalPoint>& collection_points) const;
387  void induce_signal(const PSimHit& hit,
388  const size_t hitIndex,
389  const unsigned int tofBin,
390  const PixelGeomDetUnit *pixdet,
391  const std::vector<SignalPoint>& collection_points);
392  void fluctuateEloss(int particleId, float momentum, float eloss,
393  float length, int NumberOfSegments,
394  float elossVector[],
395  CLHEP::HepRandomEngine*) const;
396  void add_noise(const PixelGeomDetUnit *pixdet,
397  float thePixelThreshold,
398  CLHEP::HepRandomEngine*);
399  void make_digis(float thePixelThresholdInE,
400  uint32_t detID,
401  const PixelGeomDetUnit* pixdet,
402  std::vector<PixelDigi>& digis,
403  std::vector<PixelDigiSimLink>& simlinks,
404  const TrackerTopology *tTopo) const;
405  void pixel_inefficiency(const PixelEfficiencies& eff,
406  const PixelGeomDetUnit* pixdet,
407  const TrackerTopology *tTopo,
408  CLHEP::HepRandomEngine*);
409 
410  void pixel_inefficiency_db(uint32_t detID);
411 
412  float pixel_aging(const PixelAging& aging,
413  const PixelGeomDetUnit* pixdet,
414  const TrackerTopology *tTopo) const;
415 
416  // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
417  const std::unique_ptr<SiPixelGainCalibrationOfflineSimService> theSiPixelGainCalibrationService_;
418  float missCalibrate(uint32_t detID, const PixelGeomDetUnit* pixdet, int col, int row, float amp) const;
419  LocalVector DriftDirection(const PixelGeomDetUnit* pixdet,
420  const GlobalVector& bfield,
421  const DetId& detId) const;
422 
423  void module_killing_conf(uint32_t detID); // remove dead modules using the list in the configuration file PixelDigi_cfi.py
424  void module_killing_DB(uint32_t detID); // remove dead modules uisng the list in the DB
425 
428 
429  double calcQ(float x) const {
430  // need erf(x/sqrt2)
431  //float x2=0.5*x*x;
432  //float a=0.147;
433  //double erf=sqrt(1.0f-exp( -1.0f*x2*( (4/M_PI)+a*x2)/(1.0+a*x2)));
434  //if (x<0.) erf*=-1.0;
435  //return 0.5*(1.0-erf);
436 
437  auto xx=std::min(0.5f*x*x,12.5f);
438  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));
439  }
440 
441 
442 };
443 
444 #endif
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)
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate
SignalPoint(float x, float y, float sigma_x, float sigma_y, float t, float a=1.0)
LocalVector drift(const StripGeomDetUnit *, const MagneticField &, const SiStripLorentzAngle &)
Definition: ShallowTools.cc:38
std::vector< std::vector< double > > thePUEfficiency
edm::ESHandle< TrackerGeometry > geom_
const std::map< int, CalParameters, std::less< int > > calmap
T sqrt(T t)
Definition: SSEVec.h:18
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_
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser
Definition: DetId.h:18
std::map< int, Amplitude, std::less< int > > signal_map_type
edm::ESHandle< SiPixelQuality > SiPixelBadModule_
GloballyPositioned< double > Frame
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:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
col
Definition: cuy.py:1008
std::map< uint32_t, signal_map_type > signalMaps
Definition: aging.py:1