CMS 3D CMS Logo

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