CMS 3D CMS Logo

Phase2TrackerDigitizerAlgorithm.h
Go to the documentation of this file.
1 #ifndef __SimTracker_SiPhase2Digitizer_Phase2TrackerDigitizerAlgorithm_h
2 #define __SimTracker_SiPhase2Digitizer_Phase2TrackerDigitizerAlgorithm_h
3 
4 #include <map>
5 #include <memory>
6 #include <vector>
7 #include <iostream>
18 
21 
22 // forward declarations
23 // For the random numbers
24 namespace CLHEP {
25  class HepRandomEngine;
26  class RandGaussQ;
27  class RandFlat;
28 }
29 
30 namespace edm {
31  class EventSetup;
32  class ParameterSet;
33 }
34 
35 class DetId;
41 class SiPixelQuality;
42 class TrackerGeometry;
43 class TrackerTopology;
44 
46  public:
47  Phase2TrackerDigitizerAlgorithm(const edm::ParameterSet& conf_common, const edm::ParameterSet& conf_specific, CLHEP::HepRandomEngine&);
49 
50  // initialization that cannot be done in the constructor
51  virtual void init(const edm::EventSetup& es) = 0;
52  virtual void initializeEvent() {
53  _signal.clear();
54  }
55  // run the algorithm to digitize a single det
56  virtual void accumulateSimHits(const std::vector<PSimHit>::const_iterator inputBegin,
57  const std::vector<PSimHit>::const_iterator inputEnd,
58  const size_t inputBeginGlobalIndex,
59  const unsigned int tofBin,
60  const Phase2TrackerGeomDetUnit* pixdet,
61  const GlobalVector& bfield) = 0;
62  virtual void digitize(const Phase2TrackerGeomDetUnit* pixdet,
63  std::map<int, DigitizerUtility::DigiSimInfo>& digi_map,
64  const TrackerTopology* tTopo);
65 
66  protected:
67  // Accessing Lorentz angle from DB:
69 
70  // Accessing Dead pixel modules from DB:
72 
73  // Accessing Map and Geom:
78  std::vector<double> barrel_efficiencies;
79  std::vector<double> endcap_efficiencies;
80  };
81 
82  // Internal type aliases
83  using signal_map_type = std::map<int, DigitizerUtility::Amplitude, std::less<int> >; // from Digi.Skel.
84  using signal_map_iterator = signal_map_type::iterator; // from Digi.Skel.
85  using signal_map_const_iterator = signal_map_type::const_iterator; // from Digi.Skel.
86  using simlink_map = std::map<unsigned int, std::vector<float>,std::less<unsigned int> > ;
87  using signalMaps = std::map<uint32_t, signal_map_type> ;
89  using Parameters = std::vector<edm::ParameterSet> ;
90 
91  // Contains the accumulated hit info.
93 
94  const bool makeDigiSimLinks_;
95 
96  const bool use_ineff_from_db_;
97  const bool use_module_killing_; // remove or not the dead pixel modules
98  const bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
99  const bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
100 
102 
103  // Variables
104  // external parameters
105  // go from Geant energy GeV to number of electrons
106  const float GeVperElectron; // 3.7E-09
107 
108  //-- drift
109  const bool alpha2Order; // Switch on/off of E.B effect
110  const bool addXtalk;
111  const float interstripCoupling;
112  const float Sigma0; //=0.0007 // Charge diffusion in microns for 300 micron Si
113  const float SigmaCoeff; // delta in the diffusion across the strip pitch
114 
115  //-- induce_signal
116  const float ClusterWidth; // Gaussian charge cutoff width in sigma units
117 
118  //-- make_digis
119  const bool doDigitalReadout; // Flag to decide analog or digital readout
120  const float theElectronPerADC; // Gain, number of electrons per adc count.
121  const int theAdcFullScale; // Saturation count, 255=8bit.
122  const float theNoiseInElectrons; // Noise (RMS) in units of electrons.
123  const float theReadoutNoise; // Noise of the readount chain in elec,
124 
125  // inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
126  const float theThresholdInE_Endcap; // threshold in electrons Endcap.
127  const float theThresholdInE_Barrel; // threshold in electrons Barrel.
128 
131 
134 
135  const float theTofLowerCut; // Cut on the particle TOF
136  const float theTofUpperCut; // Cut on the particle TOF
137  const float tanLorentzAnglePerTesla_Endcap; //FPix Lorentz angle tangent per Tesla
138  const float tanLorentzAnglePerTesla_Barrel; //BPix Lorentz angle tangent per Tesla
139 
140  // -- add_noise
141  const bool addNoise;
142  const bool addNoisyPixels;
143  const bool fluctuateCharge;
144 
145  //-- pixel efficiency
146  const bool AddPixelInefficiency; // bool to read in inefficiencies
147 
149 
150  // pseudoRadDamage
151  const double pseudoRadDamage; // Decrease the amount off freed charge that reaches the collector
152  const double pseudoRadDamageRadius; // Only apply pseudoRadDamage to pixels with radius<=pseudoRadDamageRadius
153 
154  // The PDTable
155  // HepPDTable *particleTable;
156  // ParticleDataTable *particleTable;
157 
158  //-- charge fluctuation
159  const double tMax; // The delta production cut, should be as in OSCAR = 30keV
160 
161  // Bad Pixels to be killed
162  std::vector<edm::ParameterSet> badPixels;
163 
164  // The eloss fluctuation class from G4. Is the right place?
165  const std::unique_ptr<SiG4UniversalFluctuation> fluctuate; // make a pointer
166  const std::unique_ptr<GaussianTailNoiseGenerator> theNoiser;
167 
168  //-- additional member functions
169  // Private methods
170  void primary_ionization( const PSimHit& hit, std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points) const;
171  void drift(const PSimHit& hit,
172  const Phase2TrackerGeomDetUnit* pixdet,
173  const GlobalVector& bfield,
174  const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points,
175  std::vector<DigitizerUtility::SignalPoint>& collection_points) const;
176  void induce_signal(const PSimHit& hit,
177  const size_t hitIndex,
178  const unsigned int tofBin,
179  const Phase2TrackerGeomDetUnit* pixdet,
180  const std::vector<DigitizerUtility::SignalPoint>& collection_points);
181  void fluctuateEloss(int particleId, float momentum, float eloss,
182  float length, int NumberOfSegments,
183  std::vector<float> & elossVector) const;
184  virtual void add_noise(const Phase2TrackerGeomDetUnit* pixdet, float thePixelThreshold);
185  virtual void pixel_inefficiency(const SubdetEfficiencies& eff,
186  const Phase2TrackerGeomDetUnit* pixdet,
187  const TrackerTopology* tTopo);
188 
189  virtual void pixel_inefficiency_db(uint32_t detID);
190 
191  // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
192  const std::unique_ptr<SiPixelGainCalibrationOfflineSimService> theSiPixelGainCalibrationService_;
193  LocalVector DriftDirection(const Phase2TrackerGeomDetUnit* pixdet,
194  const GlobalVector& bfield,
195  const DetId& detId) const;
196 
197  virtual void module_killing_conf(uint32_t detID); // remove dead modules using the list in the configuration file PixelDigi_cfi.py
198  virtual void module_killing_DB(uint32_t detID); // remove dead modules uisng the list in the DB
199 
201 
202  // For random numbers
203  const std::unique_ptr<CLHEP::RandFlat> flatDistribution_;
204  const std::unique_ptr<CLHEP::RandGaussQ> gaussDistribution_;
205 
206  // Threshold gaussian smearing:
207  const std::unique_ptr<CLHEP::RandGaussQ> smearedThreshold_Endcap_;
208  const std::unique_ptr<CLHEP::RandGaussQ> smearedThreshold_Barrel_;
209 
210  //for engine passed into the constructor from Digitizer
211  CLHEP::HepRandomEngine* rengine_;
212 
213  double calcQ(float x) const {
214  auto xx = std::min(0.5f * x * x,12.5f);
215  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));
216  }
217  bool pixelFlag;
218 };
219 #endif
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate
int init
Definition: HydjetWrapper.h:67
LocalVector drift(const StripGeomDetUnit *, const MagneticField &, const SiStripLorentzAngle &)
Definition: ShallowTools.cc:38
edm::ESHandle< SiPixelFedCablingMap > map_
const std::unique_ptr< SiPixelGainCalibrationOfflineSimService > theSiPixelGainCalibrationService_
std::map< unsigned int, std::vector< float >, std::less< unsigned int > > simlink_map
std::vector< edm::ParameterSet > badPixels
T sqrt(T t)
Definition: SSEVec.h:18
edm::ESHandle< SiPixelLorentzAngle > SiPixelLorentzAngle_
double f[11][100]
const std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Barrel_
T min(T a, T b)
Definition: MathUtil.h:58
const std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Endcap_
std::map< int, DigitizerUtility::Amplitude, std::less< int > > signal_map_type
signal_map_type::const_iterator signal_map_const_iterator
std::vector< edm::ParameterSet > Parameters
edm::ESHandle< TrackerGeometry > geom_
Definition: DetId.h:18
HLT enums.
const std::unique_ptr< CLHEP::RandFlat > flatDistribution_
const std::unique_ptr< CLHEP::RandGaussQ > gaussDistribution_
std::map< uint32_t, signal_map_type > signalMaps
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser
edm::ESHandle< SiPixelQuality > SiPixelBadModule_