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;
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*);
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), _hitInfo() {
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), _hitInfo(new SimHitInfoForLinks(hitp, hitIndex, tofBin) ) {
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  _hitInfo->trackIds_.pop_back();
112  }
113  }
114 
115  // can be used as a float by convers.
116  operator float() const { return _amp;}
117  float ampl() const {return _amp;}
118  std::vector<float> individualampl() const {return _frac;}
119  const std::vector<unsigned int>& trackIds() const {
120  return _hitInfo->trackIds_;
121  }
122  const std::shared_ptr<SimHitInfoForLinks>& hitInfo() const {return _hitInfo;}
123 
124  void operator+=( const Amplitude& other) {
125  _amp += other._amp;
126  //in case of contribution of noise to the digi
127  //the MC information are removed
128  if (other._frac[0]>-0.5){
129  if(other._hitInfo) {
130  std::vector<unsigned int>& otherTrackIds = other._hitInfo->trackIds_;
131  if(_hitInfo) {
132  std::vector<unsigned int>& trackIds = _hitInfo->trackIds_;
133  trackIds.insert(trackIds.end(), otherTrackIds.begin(), otherTrackIds.end());
134  } else {
135  _hitInfo.reset(new SimHitInfoForLinks(*other._hitInfo));
136  }
137  }
138  _frac.insert(_frac.end(), other._frac.begin(), other._frac.end());
139  }
140  }
141  const EncodedEventId& eventId() const {
142  return _hitInfo->eventId_;
143  }
144  const unsigned int hitIndex() const {
145  return _hitInfo->hitIndex_;
146  }
147  const unsigned int tofBin() const {
148  return _hitInfo->tofBin_;
149  }
150  void operator+=( const float& amp) {
151  _amp += amp;
152  }
153 
154  void set (const float amplitude) { // Used to reset the amplitude
155  _amp = amplitude;
156  }
157 /* void setind (const float indamplitude) { // Used to reset the amplitude */
158 /* _frac = idamplitude; */
159 /* } */
160  private:
161  float _amp;
162  std::vector<float> _frac;
163  std::shared_ptr<SimHitInfoForLinks> _hitInfo;
164  }; // end class Amplitude
165 
166  // Define a class to hold the calibration parameters per pixel
167  // Internal
169  public:
170  float p0;
171  float p1;
172  float p2;
173  float p3;
174  };
175  //
176  // Define a class for 3D ionization points and energy
177  //
182  public:
184  EnergyDepositUnit(float energy,float x, float y, float z):
185  _energy(energy),_position(x,y,z){}
187  _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  private:
193  float _energy;
195  };
196 
197  //
198  // define class to store signals on the collection surface
199  //
204  class SignalPoint {
205  public:
206  SignalPoint() : _pos(0,0), _time(0), _amplitude(0),
207  _sigma_x(1.), _sigma_y(1.), _hitp(0) {}
208 
209  SignalPoint( float x, float y, float sigma_x, float sigma_y,
210  float t, float a=1.0) :
211  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
212  _sigma_y(sigma_y), _hitp(0) {}
213 
214  SignalPoint( float x, float y, float sigma_x, float sigma_y,
215  float t, const PSimHit& hit, float a=1.0) :
216  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
217  _sigma_y(sigma_y),_hitp(&hit) {}
218 
219  const LocalPoint& position() const { return _pos;}
220  float x() const { return _pos.x();}
221  float y() const { return _pos.y();}
222  float sigma_x() const { return _sigma_x;}
223  float sigma_y() const { return _sigma_y;}
224  float time() const { return _time;}
225  float amplitude() const { return _amplitude;}
226  const PSimHit& hit() { return *_hitp;}
227  SignalPoint& set_amplitude( float amp) { _amplitude = amp; return *this;}
228 
229 
230 
231  private:
233  float _time;
234  float _amplitude;
235  float _sigma_x; // gaussian sigma in the x direction (cm)
236  float _sigma_y; // " " y direction (cm) */
237  const PSimHit* _hitp;
238  };
239 
240  //
241  // PixelEfficiencies struct
242  //
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, double> ColGeomFactors;
267  std::map<uint32_t, double> ChipGeomFactors;
268  std::map<uint32_t, size_t > iPU;
269 
271  bool matches(const DetId&, const DetId&, const std::vector<uint32_t >&);
272  };
273 
274  //
275  // PixelAging struct
276  //
280  struct PixelAging {
282  float thePixelPseudoRadDamage[20]; // PseudoRadiation Damage Values for aging studies
283  unsigned int FPixIndex; // The Efficiency index for FPix Disks
284  };
285 
286  private:
287  // Internal typedefs
288  typedef std::map<int, Amplitude, std::less<int> > signal_map_type; // from Digi.Skel.
289  typedef signal_map_type::iterator signal_map_iterator; // from Digi.Skel.
290  typedef signal_map_type::const_iterator signal_map_const_iterator; // from Digi.Skel.
291  typedef std::map<unsigned int, std::vector<float>,std::less<unsigned int> > simlink_map;
292  typedef std::map<uint32_t, signal_map_type> signalMaps;
294  typedef std::vector<edm::ParameterSet> Parameters;
295 
296  // Contains the accumulated hit info.
298 
299  const bool makeDigiSimLinks_;
300 
301  const bool use_ineff_from_db_;
302  const bool use_module_killing_; // remove or not the dead pixel modules
303  const bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
304  const bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
305 
307 
308  // Variables
309  //external parameters
310  // go from Geant energy GeV to number of electrons
311  const float GeVperElectron; // 3.7E-09
312 
313  //-- drift
314  const float Sigma0; //=0.0007 // Charge diffusion in microns for 300 micron Si
315  const float Dist300; //=0.0300 // Define 300microns for normalization
316  const bool alpha2Order; // Switch on/off of E.B effect
317 
318 
319 
320  //-- induce_signal
321  const float ClusterWidth; // Gaussian charge cutoff width in sigma units
322  //-- Allow for upgrades
323  const int NumberOfBarrelLayers; // Default = 3
324  const int NumberOfEndcapDisks; // Default = 2
325 
326  //-- make_digis
327  const float theElectronPerADC; // Gain, number of electrons per adc count.
328  const int theAdcFullScale; // Saturation count, 255=8bit.
329  const int theAdcFullScaleStack; // Saturation count for stack layers, 1=1bit.
330  const float theNoiseInElectrons; // Noise (RMS) in units of electrons.
331  const float theReadoutNoise; // Noise of the readount chain in elec,
332  //inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
333 
334  const float theThresholdInE_FPix; // Pixel threshold in electrons FPix.
335  const float theThresholdInE_BPix; // Pixel threshold in electrons BPix.
336  const float theThresholdInE_BPix_L1; // In case the BPix layer1 gets a different threshold
337 
341 
342  const double electronsPerVCAL; // for electrons - VCAL conversion
343  const double electronsPerVCAL_Offset; // in misscalibrate()
344 
345  const float theTofLowerCut; // Cut on the particle TOF
346  const float theTofUpperCut; // Cut on the particle TOF
347  const float tanLorentzAnglePerTesla_FPix; //FPix Lorentz angle tangent per Tesla
348  const float tanLorentzAnglePerTesla_BPix; //BPix Lorentz angle tangent per Tesla
349 
350  const float FPix_p0;
351  const float FPix_p1;
352  const float FPix_p2;
353  const float FPix_p3;
354  const float BPix_p0;
355  const float BPix_p1;
356  const float BPix_p2;
357  const float BPix_p3;
358 
359 
360  //-- add_noise
361  const bool addNoise;
363  const bool addNoisyPixels;
364  const bool fluctuateCharge;
365  //-- pixel efficiency
366  const bool AddPixelInefficiency; // bool to read in inefficiencies
367 
369 
370  //-- calibration smearing
371  const bool doMissCalibrate; // Switch on the calibration smearing
372  const float theGainSmearing; // The sigma of the gain fluctuation (around 1)
373  const float theOffsetSmearing; // The sigma of the offset fluct. (around 0)
374 
375  // pixel aging
376  const bool AddPixelAging;
377 
378  // The PDTable
379  //HepPDTable *particleTable;
380  //ParticleDataTable *particleTable;
381 
382  //-- charge fluctuation
383  const double tMax; // The delta production cut, should be as in OSCAR = 30keV
384  // cmsim = 100keV
385 
386  // The eloss fluctuation class from G4. Is the right place?
387  const std::unique_ptr<SiG4UniversalFluctuation> fluctuate; // make a pointer
388  const std::unique_ptr<GaussianTailNoiseGenerator> theNoiser;
389 
390  // To store calibration constants
391  const std::map<int,CalParameters,std::less<int> > calmap;
392 
393 
394  //-- additional member functions
395  // Private methods
396  std::map<int,CalParameters,std::less<int> > initCal() const;
397  void primary_ionization( const PSimHit& hit, std::vector<EnergyDepositUnit>& ionization_points, CLHEP::HepRandomEngine*) const;
398  void drift(const PSimHit& hit,
399  const PixelGeomDetUnit *pixdet,
400  const GlobalVector& bfield,
401  const TrackerTopology *tTopo,
402  const std::vector<EnergyDepositUnit>& ionization_points,
403  std::vector<SignalPoint>& collection_points) const;
404  void induce_signal(const PSimHit& hit,
405  const size_t hitIndex,
406  const unsigned int tofBin,
407  const PixelGeomDetUnit *pixdet,
408  const std::vector<SignalPoint>& collection_points);
409  void fluctuateEloss(int particleId, float momentum, float eloss,
410  float length, int NumberOfSegments,
411  float elossVector[],
412  CLHEP::HepRandomEngine*) const;
413  void add_noise(const PixelGeomDetUnit *pixdet,
414  float thePixelThreshold,
415  CLHEP::HepRandomEngine*);
416  void make_digis(float thePixelThresholdInE,
417  uint32_t detID,
418  const PixelGeomDetUnit* pixdet,
419  std::vector<PixelDigi>& digis,
420  std::vector<PixelDigiSimLink>& simlinks,
421  const TrackerTopology *tTopo) const;
422  void pixel_inefficiency(const PixelEfficiencies& eff,
423  const PixelGeomDetUnit* pixdet,
424  const TrackerTopology *tTopo,
425  CLHEP::HepRandomEngine*);
426 
427  void pixel_inefficiency_db(uint32_t detID);
428 
429  float pixel_aging(const PixelAging& aging,
430  const PixelGeomDetUnit* pixdet,
431  const TrackerTopology *tTopo) const;
432 
433  // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
434  const std::unique_ptr<SiPixelGainCalibrationOfflineSimService> theSiPixelGainCalibrationService_;
435  float missCalibrate(uint32_t detID, const PixelGeomDetUnit* pixdet, int col, int row, float amp) const;
437  const GlobalVector& bfield,
438  const DetId& detId) const;
439 
440  void module_killing_conf(uint32_t detID); // remove dead modules using the list in the configuration file PixelDigi_cfi.py
441  void module_killing_DB(uint32_t detID); // remove dead modules uisng the list in the DB
442 
445 
446  double calcQ(float x) const {
447  // need erf(x/sqrt2)
448  //float x2=0.5*x*x;
449  //float a=0.147;
450  //double erf=sqrt(1.0f-exp( -1.0f*x2*( (4/M_PI)+a*x2)/(1.0+a*x2)));
451  //if (x<0.) erf*=-1.0;
452  //return 0.5*(1.0-erf);
453 
454  auto xx=std::min(0.5f*x*x,12.5f);
455  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));
456  }
457 
458 
459 };
460 
461 #endif
void init(const edm::EventSetup &es)
Amplitude(float amp, const PSimHit *hitp, size_t hitIndex, unsigned int tofBin, float frac)
tuple t
Definition: tree.py:139
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_
void init_DynIneffDB(const edm::EventSetup &, const unsigned int &)
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
std::vector< std::vector< double > > thePUEfficiency
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
edm::ESHandle< SiPixelQuality > SiPixelBadModule_
void primary_ionization(const PSimHit &hit, std::vector< EnergyDepositUnit > &ionization_points, CLHEP::HepRandomEngine *) const
void init_from_db(const edm::ESHandle< TrackerGeometry > &, const edm::ESHandle< SiPixelDynamicInefficiency > &)
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
edm::ESHandle< SiPixelDynamicInefficiency > SiPixelDynamicInefficiency_
bool matches(const DetId &, const DetId &, const std::vector< uint32_t > &)
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)
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 *)