test
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 PixelGeomDetUnit *pixdet,
59  const GlobalVector& bfield,
60  CLHEP::HepRandomEngine*);
61  void digitize(const PixelGeomDetUnit *pixdet,
62  std::vector<PixelDigi>& digis,
63  std::vector<PixelDigiSimLink>& simlinks,
64  const TrackerTopology *tTopo,
65  CLHEP::HepRandomEngine*);
67 
68  private:
69 
70  //Accessing Lorentz angle from DB:
72 
73  //Accessing Dead pixel modules from DB:
75 
76  //Accessing Map and Geom:
79 
80  // Define internal classes
81 
82  // definition class
83  //
84  class Amplitude {
85  public:
86  Amplitude() : _amp(0.0) {}
87  Amplitude( float amp, float frac) :
88  _amp(amp), _frac(1, frac), _hitInfo() {
89  //in case of digi from noisypixels
90  //the MC information are removed
91  if (_frac[0]<-0.5) {
92  _frac.pop_back();
93  }
94  }
95 
96  Amplitude( float amp, const PSimHit* hitp, float frac) :
97  _amp(amp), _frac(1, frac), _hitInfo(new SimHitInfoForLinks(hitp)) {
98 
99  //in case of digi from noisypixels
100  //the MC information are removed
101  if (_frac[0]<-0.5) {
102  _frac.pop_back();
103  _hitInfo->trackIds_.pop_back();
104  }
105  }
106 
107  // can be used as a float by convers.
108  operator float() const { return _amp;}
109  float ampl() const {return _amp;}
110  std::vector<float> individualampl() const {return _frac;}
111  const std::vector<unsigned int>& trackIds() const {
112  return _hitInfo->trackIds_;
113  }
114  const std::shared_ptr<SimHitInfoForLinks>& hitInfo() const {return _hitInfo;}
115 
116  void operator+=( const Amplitude& other) {
117  _amp += other._amp;
118  //in case of contribution of noise to the digi
119  //the MC information are removed
120  if (other._frac[0]>-0.5){
121  if(other._hitInfo) {
122  std::vector<unsigned int>& otherTrackIds = other._hitInfo->trackIds_;
123  if(_hitInfo) {
124  std::vector<unsigned int>& trackIds = _hitInfo->trackIds_;
125  trackIds.insert(trackIds.end(), otherTrackIds.begin(), otherTrackIds.end());
126  } else {
127  _hitInfo.reset(new SimHitInfoForLinks(*other._hitInfo));
128  }
129  }
130  _frac.insert(_frac.end(), other._frac.begin(), other._frac.end());
131  }
132  }
133  const EncodedEventId& eventId() const {
134  return _hitInfo->eventId_;
135  }
136  void operator+=( const float& amp) {
137  _amp += amp;
138  }
139 
140  void set (const float amplitude) { // Used to reset the amplitude
141  _amp = amplitude;
142  }
143 /* void setind (const float indamplitude) { // Used to reset the amplitude */
144 /* _frac = idamplitude; */
145 /* } */
146  private:
147  float _amp;
148  std::vector<float> _frac;
149  std::shared_ptr<SimHitInfoForLinks> _hitInfo;
150  }; // end class Amplitude
151 
152  // Define a class to hold the calibration parameters per pixel
153  // Internal
155  public:
156  float p0;
157  float p1;
158  float p2;
159  float p3;
160  };
161  //
162  // Define a class for 3D ionization points and energy
163  //
168  public:
170  EnergyDepositUnit(float energy,float x, float y, float z):
171  _energy(energy),_position(x,y,z){}
173  _energy(energy),_position(position){}
174  float x() const{return _position.x();}
175  float y() const{return _position.y();}
176  float z() const{return _position.z();}
177  float energy() const { return _energy;}
178  private:
179  float _energy;
181  };
182 
183  //
184  // define class to store signals on the collection surface
185  //
190  class SignalPoint {
191  public:
192  SignalPoint() : _pos(0,0), _time(0), _amplitude(0),
193  _sigma_x(1.), _sigma_y(1.), _hitp(0) {}
194 
195  SignalPoint( float x, float y, float sigma_x, float sigma_y,
196  float t, float a=1.0) :
197  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
198  _sigma_y(sigma_y), _hitp(0) {}
199 
200  SignalPoint( float x, float y, float sigma_x, float sigma_y,
201  float t, const PSimHit& hit, float a=1.0) :
202  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
203  _sigma_y(sigma_y),_hitp(&hit) {}
204 
205  const LocalPoint& position() const { return _pos;}
206  float x() const { return _pos.x();}
207  float y() const { return _pos.y();}
208  float sigma_x() const { return _sigma_x;}
209  float sigma_y() const { return _sigma_y;}
210  float time() const { return _time;}
211  float amplitude() const { return _amplitude;}
212  const PSimHit& hit() { return *_hitp;}
213  SignalPoint& set_amplitude( float amp) { _amplitude = amp; return *this;}
214 
215 
216 
217  private:
219  float _time;
220  float _amplitude;
221  float _sigma_x; // gaussian sigma in the x direction (cm)
222  float _sigma_y; // " " y direction (cm) */
223  const PSimHit* _hitp;
224  };
225 
226  //
227  // PixelEfficiencies struct
228  //
234  float thePixelEfficiency[20]; // Single pixel effciency
235  float thePixelColEfficiency[20]; // Column effciency
236  float thePixelChipEfficiency[20]; // ROC efficiency
237  std::vector<double> theLadderEfficiency_BPix[20]; // Ladder efficiency
238  std::vector<double> theModuleEfficiency_BPix[20]; // Module efficiency
239  std::vector<double> thePUEfficiency_BPix[20]; // Instlumi dependent efficiency
240  unsigned int FPixIndex; // The Efficiency index for FPix Disks
241  };
242 
243  private:
244  // Needed by dynamic inefficiency
245  // 0-3 BPix, 4-5 FPix
246  double _pu_scale[20];
247 
248  // Internal typedefs
249  typedef std::map<int, Amplitude, std::less<int> > signal_map_type; // from Digi.Skel.
250  typedef signal_map_type::iterator signal_map_iterator; // from Digi.Skel.
251  typedef signal_map_type::const_iterator signal_map_const_iterator; // from Digi.Skel.
252  typedef std::map<unsigned int, std::vector<float>,std::less<unsigned int> > simlink_map;
253  typedef std::map<uint32_t, signal_map_type> signalMaps;
255  typedef std::vector<edm::ParameterSet> Parameters;
256 
257  // Contains the accumulated hit info.
259 
260  const bool makeDigiSimLinks_;
261 
262  const bool use_ineff_from_db_;
263  const bool use_module_killing_; // remove or not the dead pixel modules
264  const bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
265  const bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
266 
268 
269  // Variables
270  //external parameters
271  // go from Geant energy GeV to number of electrons
272  const float GeVperElectron; // 3.7E-09
273 
274  //-- drift
275  const float Sigma0; //=0.0007 // Charge diffusion in microns for 300 micron Si
276  const float Dist300; //=0.0300 // Define 300microns for normalization
277  const bool alpha2Order; // Switch on/off of E.B effect
278 
279 
280 
281  //-- induce_signal
282  const float ClusterWidth; // Gaussian charge cutoff width in sigma units
283  //-- Allow for upgrades
284  const int NumberOfBarrelLayers; // Default = 3
285  const int NumberOfEndcapDisks; // Default = 2
286 
287  //-- make_digis
288  const float theElectronPerADC; // Gain, number of electrons per adc count.
289  const int theAdcFullScale; // Saturation count, 255=8bit.
290  const int theAdcFullScaleStack; // Saturation count for stack layers, 1=1bit.
291  const int theFirstStackLayer; // The first BPix layer to use theAdcFullScaleStack.
292  const float theNoiseInElectrons; // Noise (RMS) in units of electrons.
293  const float theReadoutNoise; // Noise of the readount chain in elec,
294  //inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
295 
296  const float theThresholdInE_FPix; // Pixel threshold in electrons FPix.
297  const float theThresholdInE_BPix; // Pixel threshold in electrons BPix.
298  const float theThresholdInE_BPix_L1; // In case the BPix layer1 gets a different threshold
299 
303 
304  const double electronsPerVCAL; // for electrons - VCAL conversion
305  const double electronsPerVCAL_Offset; // in misscalibrate()
306 
307  const float theTofLowerCut; // Cut on the particle TOF
308  const float theTofUpperCut; // Cut on the particle TOF
309  const float tanLorentzAnglePerTesla_FPix; //FPix Lorentz angle tangent per Tesla
310  const float tanLorentzAnglePerTesla_BPix; //BPix Lorentz angle tangent per Tesla
311 
312  const float FPix_p0;
313  const float FPix_p1;
314  const float FPix_p2;
315  const float FPix_p3;
316  const float BPix_p0;
317  const float BPix_p1;
318  const float BPix_p2;
319  const float BPix_p3;
320 
321 
322  //-- add_noise
323  const bool addNoise;
325  const bool addNoisyPixels;
326  const bool fluctuateCharge;
327  //-- pixel efficiency
328  const bool AddPixelInefficiency; // bool to read in inefficiencies
329 
331 
332  //-- calibration smearing
333  const bool doMissCalibrate; // Switch on the calibration smearing
334  const float theGainSmearing; // The sigma of the gain fluctuation (around 1)
335  const float theOffsetSmearing; // The sigma of the offset fluct. (around 0)
336 
337  // pseudoRadDamage
338  const double pseudoRadDamage; // Decrease the amount off freed charge that reaches the collector
339  const double pseudoRadDamageRadius; // Only apply pseudoRadDamage to pixels with radius<=pseudoRadDamageRadius
340  // The PDTable
341  //HepPDTable *particleTable;
342  //ParticleDataTable *particleTable;
343 
344  //-- charge fluctuation
345  const double tMax; // The delta production cut, should be as in OSCAR = 30keV
346  // cmsim = 100keV
347 
348  // The eloss fluctuation class from G4. Is the right place?
349  const std::unique_ptr<SiG4UniversalFluctuation> fluctuate; // make a pointer
350  const std::unique_ptr<GaussianTailNoiseGenerator> theNoiser;
351 
352  // To store calibration constants
353  const std::map<int,CalParameters,std::less<int> > calmap;
354 
355 
356  //-- additional member functions
357  // Private methods
358  std::map<int,CalParameters,std::less<int> > initCal() const;
359  void primary_ionization( const PSimHit& hit, std::vector<EnergyDepositUnit>& ionization_points, CLHEP::HepRandomEngine*) const;
360  void drift(const PSimHit& hit,
361  const PixelGeomDetUnit *pixdet,
362  const GlobalVector& bfield,
363  const std::vector<EnergyDepositUnit>& ionization_points,
364  std::vector<SignalPoint>& collection_points) const;
365  void induce_signal(const PSimHit& hit,
366  const PixelGeomDetUnit *pixdet,
367  const std::vector<SignalPoint>& collection_points);
368  void fluctuateEloss(int particleId, float momentum, float eloss,
369  float length, int NumberOfSegments,
370  float elossVector[],
371  CLHEP::HepRandomEngine*) const;
372  void add_noise(const PixelGeomDetUnit *pixdet,
373  float thePixelThreshold,
374  CLHEP::HepRandomEngine*);
375  void make_digis(float thePixelThresholdInE,
376  uint32_t detID,
377  std::vector<PixelDigi>& digis,
378  std::vector<PixelDigiSimLink>& simlinks,
379  const TrackerTopology *tTopo) const;
380  void pixel_inefficiency(const PixelEfficiencies& eff,
381  const PixelGeomDetUnit* pixdet,
382  const TrackerTopology *tTopo,
383  CLHEP::HepRandomEngine*);
384 
385  void pixel_inefficiency_db(uint32_t detID);
386 
387  // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
388  const std::unique_ptr<SiPixelGainCalibrationOfflineSimService> theSiPixelGainCalibrationService_;
389  float missCalibrate(uint32_t detID, int col, int row, float amp) const;
391  const GlobalVector& bfield,
392  const DetId& detId) const;
393 
394  void module_killing_conf(uint32_t detID); // remove dead modules using the list in the configuration file PixelDigi_cfi.py
395  void module_killing_DB(uint32_t detID); // remove dead modules uisng the list in the DB
396 
398 
399  double calcQ(float x) const {
400  // need erf(x/sqrt2)
401  //float x2=0.5*x*x;
402  //float a=0.147;
403  //double erf=sqrt(1.0f-exp( -1.0f*x2*( (4/M_PI)+a*x2)/(1.0+a*x2)));
404  //if (x<0.) erf*=-1.0;
405  //return 0.5*(1.0-erf);
406 
407  auto xx=std::min(0.5f*x*x,12.5f);
408  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));
409  }
410 
411 
412 };
413 
414 #endif
void init(const edm::EventSetup &es)
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)
float missCalibrate(uint32_t detID, int col, int row, float amp) const
T y() const
Definition: PV3DBase.h:63
edm::ESHandle< TrackerGeometry > geom_
const std::shared_ptr< SimHitInfoForLinks > & hitInfo() const
Amplitude(float amp, const PSimHit *hitp, float frac)
void induce_signal(const PSimHit &hit, const PixelGeomDetUnit *pixdet, const std::vector< SignalPoint > &collection_points)
const std::map< int, CalParameters, std::less< int > > calmap
void drift(const PSimHit &hit, const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const std::vector< EnergyDepositUnit > &ionization_points, std::vector< SignalPoint > &collection_points) 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]
void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, CLHEP::HepRandomEngine *)
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
void make_digis(float thePixelThresholdInE, uint32_t detID, std::vector< PixelDigi > &digis, std::vector< PixelDigiSimLink > &simlinks, const TrackerTopology *tTopo) const
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)
GloballyPositioned< double > Frame
std::vector< edm::ParameterSet > Parameters
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
Definition: DDAxes.h:10
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 *)