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 
5 #include <string>
6 #include <map>
7 
13 
14 //#include "SimGeneral/HepPDT/interface/HepPDTable.h"
15 //#include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
16 
23 //#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeomFromDetUnits.h"
25 //#include "DataFormats/GeometrySurface/interface/TkRotation.h"
26 //#include "DataFormats/GeometrySurface/interface/GloballyPositioned.h"
28 //#include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLinkCollection.h"
30 
32 
33 // pixel gain payload access (offline version for Simulation)
35 
36 // Accessing Pixel Lorentz Angle from the DB:
39 
40 // Accessing Pixel dead modules from the DB:
43 
46 
50 
51 
52 // For the random numbers
53 namespace CLHEP {
54  class HepRandomEngine;
55  class RandGaussQ;
56  class RandFlat;
57 }
58 
60  public:
61 
62  SiPixelDigitizerAlgorithm(const edm::ParameterSet& conf, CLHEP::HepRandomEngine&);
64 
65  //run the algorithm to digitize a single det
67  run(const std::vector<PSimHit> &input,PixelGeomDetUnit *pixdet,GlobalVector);
68 
69  //
70  std::vector<PixelDigiSimLink> make_link() {
71  return link_coll; }
72  void init(const edm::EventSetup& es);
73  void fillDeadModules(const edm::EventSetup& es);
74  void fillLorentzAngle(const edm::EventSetup& es);
75  void fillMapandGeom(const edm::EventSetup& es);
76 
77 
78  private:
79 
80  //Accessing Lorentz angle from DB:
82 
83  //Accessing Dead pixel modules from DB:
85 
86  //Accessing Map and Geom:
89 
90  typedef std::vector<edm::ParameterSet> Parameters;
92 
93  // Define internal classes
94 
95  // Define a class to hold the calibration parameters per pixel
96  // Internal
97  class CalParameters {
98  public:
99  float p0;
100  float p1;
101  float p2;
102  float p3;
103  };
104  //
105  // Define a class for 3D ionization points and energy
106  //
111  public:
113  EnergyDepositUnit(float energy,float x, float y, float z):
114  _energy(energy),_position(x,y,z){}
116  _energy(energy),_position(position){}
117  float x() const{return _position.x();}
118  float y() const{return _position.y();}
119  float z() const{return _position.z();}
120  float energy() const { return _energy;}
121  private:
122  float _energy;
124  };
125 
126  //
127  // define class to store signals on the collection surface
128  //
133  class SignalPoint {
134  public:
135  SignalPoint() : _pos(0,0), _time(0), _amplitude(0),
136  _sigma_x(1.), _sigma_y(1.), _hitp(0) {}
137 
138  SignalPoint( float x, float y, float sigma_x, float sigma_y,
139  float t, float a=1.0) :
140  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
141  _sigma_y(sigma_y), _hitp(0) {}
142 
143  SignalPoint( float x, float y, float sigma_x, float sigma_y,
144  float t, const PSimHit& hit, float a=1.0) :
145  _pos(x,y), _time(t), _amplitude(a), _sigma_x(sigma_x),
146  _sigma_y(sigma_y),_hitp(&hit) {}
147 
148  const LocalPoint& position() const { return _pos;}
149  float x() const { return _pos.x();}
150  float y() const { return _pos.y();}
151  float sigma_x() const { return _sigma_x;}
152  float sigma_y() const { return _sigma_y;}
153  float time() const { return _time;}
154  float amplitude() const { return _amplitude;}
155  const PSimHit& hit() { return *_hitp;}
156  SignalPoint& set_amplitude( float amp) { _amplitude = amp; return *this;}
157 
158 
159 
160  private:
162  float _time;
163  float _amplitude;
164  float _sigma_x; // gaussian sigma in the x direction (cm)
165  float _sigma_y; // " " y direction (cm) */
166  const PSimHit* _hitp;
167  };
168 
169  //
170  // definition class
171  //
176  class Amplitude {
177  public:
178  Amplitude() : _amp(0.0) { _hits.reserve(1);}
179  Amplitude( float amp, const PSimHit* hitp, float frac) :
180  _amp(amp), _hits(1, hitp), _frac(1,frac) {
181 
182  //in case of digi from noisypixels
183  //the MC information are removed
184  if (_frac[0]<-0.5) {
185  _frac.pop_back();
186  _hits.pop_back();
187  }
188 
189  }
190 
191  // can be used as a float by convers.
192  operator float() const { return _amp;}
193  float ampl() const {return _amp;}
194  std::vector<float> individualampl() const {return _frac;}
195  const std::vector<const PSimHit*>& hits() { return _hits;}
196 
197  void operator+=( const Amplitude& other) {
198  _amp += other._amp;
199  //in case of contribution of noise to the digi
200  //the MC information are removed
201  if (other._frac[0]>-0.5){
202  _hits.insert( _hits.end(), other._hits.begin(), other._hits.end());
203  _frac.insert(_frac.end(), other._frac.begin(), other._frac.end());
204  }
205  }
206 
207  void operator+=( const float& amp) {
208  _amp += amp;
209  }
210 
211  void set (const float amplitude) { // Used to reset the amplitude
212  _amp = amplitude;
213  }
214 /* void setind (const float indamplitude) { // Used to reset the amplitude */
215 /* _frac = idamplitude; */
216 /* } */
217  private:
218  float _amp;
219  std::vector<const PSimHit*> _hits;
220  std::vector<float> _frac;
221  }; // end class Amplitude
222 
223 
224  private:
225 
226  // Internal typedefs
227  typedef std::map< int, Amplitude, std::less<int> > signal_map_type; // from
228  typedef signal_map_type::iterator signal_map_iterator; //Digi.Skel.
229  typedef std::map<unsigned int, std::vector<float>,std::less<unsigned int> >
232 
233  // Variables
235  //external parameters
236  //-- primary ionization
237  int NumberOfSegments; // =20 does not work ;
238  // go from Geant energy GeV to number of electrons
239  float GeVperElectron; // 3.7E-09
240 
241  //-- drift
242  float Sigma0; //=0.0007 // Charge diffusion in microns for 300 micron Si
243  float Dist300; //=0.0300 // Define 300microns for normalization
244  bool alpha2Order; // Switch on/off of E.B effect
245 
246 
247 
248  //-- induce_signal
249  float ClusterWidth; // Gaussian charge cutoff width in sigma units
250  //-- make_digis
251  float theElectronPerADC; // Gain, number of electrons per adc count.
252  int theAdcFullScale; // Saturation count, 255=8bit.
253  float theNoiseInElectrons; // Noise (RMS) in units of electrons.
254  float theReadoutNoise; // Noise of the readount chain in elec,
255  //inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
256 
258 
259  float thePixelThreshold; // Pixel threshold in units of noise.
260 
261  float thePixelThresholdInE; // Pixel noise in electrons.
262 
263  float theThresholdInE_FPix; // Pixel threshold in electrons FPix.
264  float theThresholdInE_BPix; // Pixel threshold in electrons BPix.
265 
268 
269  double electronsPerVCAL; // for electrons - VCAL conversion
270  double electronsPerVCAL_Offset; // in misscalibrate()
271 
272  float theTofLowerCut; // Cut on the particle TOF
273  float theTofUpperCut; // Cut on the particle TOF
274  float tanLorentzAnglePerTesla_FPix; //FPix Lorentz angle tangent per Tesla
275  float tanLorentzAnglePerTesla_BPix; //BPix Lorentz angle tangent per Tesla
276 
277  float FPix_p0;
278  float FPix_p1;
279  float FPix_p2;
280  float FPix_p3;
281  float BPix_p0;
282  float BPix_p1;
283  float BPix_p2;
284  float BPix_p3;
285 
286 
287  //-- add_noise
288  bool addNoise;
293  //-- pixel efficiency
294  bool pixelInefficiency; // Switch on pixel ineffciency
295  int thePixelLuminosity; // luminosity for inefficiency, 0,1,10
296 
298 
299  int theColsInChip; // num of columns per ROC (for pix ineff.)
300  int theRowsInChip; // num of rows per ROC
301 
302  int numColumns; // number of pixel columns in a module (detUnit)
303  int numRows; // number rows
304  float moduleThickness; // sensor thickness
305  // int digis;
307  uint32_t detID; // Det id
308 
309 
310  std::vector<PSimHit> _PixelHits; //cache
312 
313  std::vector<PixelDigi> internal_coll; //empty vector of PixelDigi used in digitize
314 
315  std::vector<PixelDigiSimLink> link_coll;
317 
318  float PixelEff;
319  float PixelColEff;
324  float thePixelEfficiency[6]; // Single pixel effciency
325  float thePixelColEfficiency[6]; // Column effciency
326  float thePixelChipEfficiency[6]; // ROC efficiency
327 
328  //-- calibration smearing
329  bool doMissCalibrate; // Switch on the calibration smearing
330  float theGainSmearing; // The sigma of the gain fluctuation (around 1)
331  float theOffsetSmearing; // The sigma of the offset fluct. (around 0)
332 
333 
334  // The PDTable
335  //HepPDTable *particleTable;
336  //ParticleDataTable *particleTable;
337 
338  //-- charge fluctuation
339  double tMax; // The delta production cut, should be as in OSCAR = 30keV
340  // cmsim = 100keV
341  // The eloss fluctuation class from G4. Is the right place?
342  SiG4UniversalFluctuation * fluctuate; // make a pointer
344 
345 
346 
347 
348  PixelIndices * pIndexConverter; // Pointer to the index converter
349 
350  std::vector<EnergyDepositUnit> _ionization_points;
351  std::vector<SignalPoint> _collection_points;
352 
354  signal_map_type _signal; // from Digi.Skel.
355 
356  // To store calibration constants
357  std::map<int,CalParameters,std::less<int> > calmap;
358 
359 
360  //-- additional member functions
361  // Private methods
362  void primary_ionization( const PSimHit& hit);
363  std::vector<PixelDigi> digitize(PixelGeomDetUnit *det);
364  void drift(const PSimHit& hit);
365  void induce_signal( const PSimHit& hit);
366  void fluctuateEloss(int particleId, float momentum, float eloss,
367  float length, int NumberOfSegments,
368  float elossVector[]);
369  void add_noise();
370  void make_digis();
371  void pixel_inefficiency();
373 
374  bool use_module_killing_; // remove or not the dead pixel modules
375  bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
376  bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
377 
378  void pixel_inefficiency_db();
379  // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
381  float missCalibrate(int col, int row, float amp) const;
383 
384  void module_killing_conf(); // remove dead modules using the list in the configuration file PixelDigi_cfi.py
385  void module_killing_DB(); // remove dead modules uisng the list in the DB
386 
387  // For random numbers
388  CLHEP::HepRandomEngine& rndEngine;
389 
390  CLHEP::RandFlat *flatDistribution_;
391  CLHEP::RandGaussQ *gaussDistribution_;
392  CLHEP::RandGaussQ *gaussDistributionVCALNoise_;
393 
394 
395  // Threshold gaussian smearing:
396  CLHEP::RandGaussQ *smearedThreshold_FPix_;
397  CLHEP::RandGaussQ *smearedThreshold_BPix_;
398 
399  CLHEP::RandGaussQ *smearedChargeDistribution_ ;
400 
401  // the random generator
402  CLHEP::RandGaussQ* theGaussianDistribution;
403 
404 
405 };
406 
407 #endif
void init(const edm::EventSetup &es)
std::vector< T > collection_type
Definition: DetSet.h:25
std::vector< float > individualampl() const
edm::ESHandle< SiPixelFedCablingMap > map_
GaussianTailNoiseGenerator * theNoiser
SignalPoint(float x, float y, float sigma_x, float sigma_y, float t, const PSimHit &hit, float a=1.0)
std::map< unsigned int, std::vector< float >, std::less< unsigned int > > simlink_map
std::vector< PixelDigiSimLink > make_link()
EnergyDepositUnit(float energy, Local3DPoint position)
const PixelGeomDetUnit * _detp
CLHEP::RandGaussQ * smearedChargeDistribution_
SignalPoint(float x, float y, float sigma_x, float sigma_y, float t, float a=1.0)
SiG4UniversalFluctuation * fluctuate
std::map< int, Amplitude, std::less< int > > signal_map_type
T y() const
Definition: PV3DBase.h:57
std::vector< PixelDigi > digitize(PixelGeomDetUnit *det)
edm::ESHandle< TrackerGeometry > geom_
void fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments, float elossVector[])
const std::vector< const PSimHit * > & hits()
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
SiPixelGainCalibrationOfflineSimService * theSiPixelGainCalibrationService_
CLHEP::RandGaussQ * gaussDistributionVCALNoise_
Amplitude(float amp, const PSimHit *hitp, float frac)
void fillDeadModules(const edm::EventSetup &es)
T z() const
Definition: PV3DBase.h:58
std::vector< PixelDigiSimLink > link_coll
edm::DetSet< PixelDigi >::collection_type run(const std::vector< PSimHit > &input, PixelGeomDetUnit *pixdet, GlobalVector)
CLHEP::RandGaussQ * theGaussianDistribution
void fillLorentzAngle(const edm::EventSetup &es)
signal_map_type::iterator signal_map_iterator
EnergyDepositUnit(float energy, float x, float y, float z)
edm::ESHandle< SiPixelLorentzAngle > SiPixelLorentzAngle_
std::map< int, CalParameters, std::less< int > > calmap
tuple conf
Definition: dbtoconf.py:185
tuple input
Definition: collect_tpl.py:10
std::vector< PixelDigi > internal_coll
edm::ESHandle< SiPixelQuality > SiPixelBadModule_
void fillMapandGeom(const edm::EventSetup &es)
GloballyPositioned< double > Frame
CLHEP::RandGaussQ * smearedThreshold_FPix_
std::vector< EnergyDepositUnit > _ionization_points
std::vector< edm::ParameterSet > Parameters
CLHEP::RandGaussQ * smearedThreshold_BPix_
double a
Definition: hdecay.h:121
void primary_ionization(const PSimHit &hit)
SiPixelDigitizerAlgorithm(const edm::ParameterSet &conf, CLHEP::HepRandomEngine &)
std::vector< SignalPoint > _collection_points
float missCalibrate(int col, int row, float amp) const
T x() const
Definition: PV3DBase.h:56
void induce_signal(const PSimHit &hit)
CLHEP::HepRandomEngine & rndEngine