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