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