CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
15 
18 
19 // Units and Constants
21 #include "CLHEP/Units/GlobalPhysicalConstants.h"
22 #include "CLHEP/Units/GlobalSystemOfUnits.h"
23 
24 // forward declarations
25 // For the random numbers
26 namespace CLHEP {
27  class HepRandomEngine;
28  class RandGaussQ;
29  class RandFlat;
30 } // namespace CLHEP
31 
32 class DetId;
38 class SiPixelQuality;
40 class TrackerGeometry;
41 class TrackerTopology;
42 
43 // REMEMBER CMS conventions:
44 // -- Energy: GeV
45 // -- momentum: GeV/c
46 // -- mass: GeV/c^2
47 // -- Distance, position: cm
48 // -- Time: ns
49 // -- Angles: radian
50 // Some constants in convenient units
51 constexpr double c_cm_ns = CLHEP::c_light * CLHEP::ns / CLHEP::cm;
52 constexpr double c_inv = 1.0 / c_cm_ns;
53 
55 public:
57  const edm::ParameterSet& conf_specific,
60 
61  // initialization that cannot be done in the constructor
62  virtual void init(const edm::EventSetup& es) = 0;
63  virtual void initializeEvent(CLHEP::HepRandomEngine& eng);
64 
65  // run the algorithm to digitize a single det
66  virtual void accumulateSimHits(const std::vector<PSimHit>::const_iterator inputBegin,
67  const std::vector<PSimHit>::const_iterator inputEnd,
68  const size_t inputBeginGlobalIndex,
69  const uint32_t tofBin,
70  const Phase2TrackerGeomDetUnit* pixdet,
71  const GlobalVector& bfield);
72  virtual void digitize(const Phase2TrackerGeomDetUnit* pixdet,
73  std::map<int, DigitizerUtility::DigiSimInfo>& digi_map,
74  const TrackerTopology* tTopo);
75  virtual bool select_hit(const PSimHit& hit, double tCorr, double& sigScale) const { return true; }
76  virtual bool isAboveThreshold(const DigitizerUtility::SimHitInfo* hitInfo, float charge, float thr) const {
77  return true;
78  }
79 
80  // For premixing
81  void loadAccumulator(uint32_t detId, const std::map<int, float>& accumulator);
82 
83 protected:
84  // Accessing Inner Tracker Lorentz angle from DB:
86 
87  // Accessing Outer Tracker Lorentz angle from DB:
89 
90  // Accessing Dead pixel modules from DB:
92 
93  // Accessing Map and Geom:
98  std::vector<double> barrel_efficiencies;
99  std::vector<double> endcap_efficiencies;
100  };
101 
102  // Internal type aliases
103  using signal_map_type = std::map<int, DigitizerUtility::Amplitude, std::less<int> >;
104  using signalMaps = std::map<uint32_t, signal_map_type>;
106  using Parameters = std::vector<edm::ParameterSet>;
107 
108  // Contains the accumulated hit info.
110 
111  const bool makeDigiSimLinks_;
112 
113  const bool use_ineff_from_db_;
114  const bool use_module_killing_; // remove or not the dead pixel modules
115  const bool use_deadmodule_DB_; // if we want to get dead pixel modules from the DataBase.
116  const bool use_LorentzAngle_DB_; // if we want to get Lorentz angle from the DataBase.
117 
119 
120  // Variables
121  // external parameters
122  // go from Geant energy GeV to number of electrons
123  const float GeVperElectron_; // 3.7E-09
124 
125  //-- drift
126  const bool alpha2Order_; // Switch on/off of E.B effect
127  const bool addXtalk_;
128  const float interstripCoupling_;
129  const float Sigma0_; //=0.0007 // Charge diffusion in microns for 300 micron Si
130  const float SigmaCoeff_; // delta in the diffusion across the strip pitch
131 
132  //-- induce_signal
133  const float clusterWidth_; // Gaussian charge cutoff width in sigma units
134 
135  //-- make_digis
136  const int thePhase2ReadoutMode_; // Flag to decide readout mode (digital/Analog dual slope etc.)
137  const float theElectronPerADC_; // Gain, number of electrons per adc count.
138  const int theAdcFullScale_; // Saturation count, 255=8bit.
139  const float theNoiseInElectrons_; // Noise (RMS) in units of electrons.
140  const float theReadoutNoise_; // Noise of the readount chain in elec,
141 
142  // inludes DCOL-Amp,TBM-Amp, Alt, AOH,OptRec.
143  const float theThresholdInE_Endcap_; // threshold in electrons Endcap.
144  const float theThresholdInE_Barrel_; // threshold in electrons Barrel.
145 
148 
151 
152  const float theTofLowerCut_; // Cut on the particle TOF
153  const float theTofUpperCut_; // Cut on the particle TOF
154  const float tanLorentzAnglePerTesla_Endcap_; //FPix Lorentz angle tangent per Tesla
155  const float tanLorentzAnglePerTesla_Barrel_; //BPix Lorentz angle tangent per Tesla
156 
157  // -- add_noise
158  const bool addNoise_;
159  const bool addNoisyPixels_;
160  const bool fluctuateCharge_;
161 
162  //-- pixel efficiency
163  const bool addPixelInefficiency_; // bool to read in inefficiencies
164 
166 
167  // pseudoRadDamage
168  const double pseudoRadDamage_; // Decrease the amount off freed charge that reaches the collector
169  const double pseudoRadDamageRadius_; // Only apply pseudoRadDamage to pixels with radius<=pseudoRadDamageRadius
170 
171  // The PDTable
172  // HepPDTable *particleTable;
173  // ParticleDataTable *particleTable;
174 
175  //-- charge fluctuation
176  const double tMax_; // The delta production cut, should be as in OSCAR = 30keV
177 
178  // Bad Pixels to be killed
180 
181  // The eloss fluctuation class from G4. Is the right place?
182  const std::unique_ptr<SiG4UniversalFluctuation> fluctuate_; // make a pointer
183  const std::unique_ptr<GaussianTailNoiseGenerator> theNoiser_;
184 
185  //-- additional member functions
186  // Private methods
187  virtual std::vector<DigitizerUtility::EnergyDepositUnit> primary_ionization(const PSimHit& hit) const;
188  virtual std::vector<DigitizerUtility::SignalPoint> drift(
189  const PSimHit& hit,
190  const Phase2TrackerGeomDetUnit* pixdet,
191  const GlobalVector& bfield,
192  const std::vector<DigitizerUtility::EnergyDepositUnit>& ionization_points) const;
193  virtual void induce_signal(const PSimHit& hit,
194  const size_t hitIndex,
195  const uint32_t tofBin,
196  const Phase2TrackerGeomDetUnit* pixdet,
197  const std::vector<DigitizerUtility::SignalPoint>& collection_points);
198  virtual std::vector<float> fluctuateEloss(
199  int particleId, float momentum, float eloss, float length, int NumberOfSegments) const;
200  virtual void add_noise(const Phase2TrackerGeomDetUnit* pixdet);
201  virtual void add_cross_talk(const Phase2TrackerGeomDetUnit* pixdet);
202  virtual void add_noisy_cells(const Phase2TrackerGeomDetUnit* pixdet, float thePixelThreshold);
203  virtual void pixel_inefficiency(const SubdetEfficiencies& eff,
204  const Phase2TrackerGeomDetUnit* pixdet,
205  const TrackerTopology* tTopo);
206 
207  virtual void pixel_inefficiency_db(uint32_t detID);
208 
209  // access to the gain calibration payloads in the db. Only gets initialized if check_dead_pixels_ is set to true.
210  const std::unique_ptr<SiPixelGainCalibrationOfflineSimService> theSiPixelGainCalibrationService_;
212  const GlobalVector& bfield,
213  const DetId& detId) const;
214 
215  // remove dead modules using the list in the configuration file PixelDigi_cfi.py
216  virtual void module_killing_conf(uint32_t detID);
217  virtual void module_killing_DB(
218  const Phase2TrackerGeomDetUnit* pixdet); // remove dead modules uisng the list in the DB
219 
221  float calcQ(float x);
222 
223  // For random numbers
224  std::unique_ptr<CLHEP::RandGaussQ> gaussDistribution_;
225 
226  // Threshold gaussian smearing:
227  std::unique_ptr<CLHEP::RandGaussQ> smearedThreshold_Endcap_;
228  std::unique_ptr<CLHEP::RandGaussQ> smearedThreshold_Barrel_;
229 
230  //for engine passed into the constructor from Digitizer
231  CLHEP::HepRandomEngine* rengine_;
232 
233  // convert signal in electrons to ADC counts
234  int convertSignalToAdc(uint32_t detID, float signal_in_elec, float threshold);
235 
237 };
238 #endif
virtual void module_killing_DB(const Phase2TrackerGeomDetUnit *pixdet)
LocalVector DriftDirection(const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const DetId &detId) const
void loadAccumulator(uint32_t detId, const std::map< int, float > &accumulator)
virtual bool select_hit(const PSimHit &hit, double tCorr, double &sigScale) const
constexpr double c_cm_ns
virtual bool isAboveThreshold(const DigitizerUtility::SimHitInfo *hitInfo, float charge, float thr) const
virtual void initializeEvent(CLHEP::HepRandomEngine &eng)
virtual void add_cross_talk(const Phase2TrackerGeomDetUnit *pixdet)
const SiPixelLorentzAngle * siPixelLorentzAngle_
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate_
virtual void add_noise(const Phase2TrackerGeomDetUnit *pixdet)
virtual void module_killing_conf(uint32_t detID)
const std::unique_ptr< SiPixelGainCalibrationOfflineSimService > theSiPixelGainCalibrationService_
virtual void add_noisy_cells(const Phase2TrackerGeomDetUnit *pixdet, float thePixelThreshold)
std::unique_ptr< CLHEP::RandGaussQ > gaussDistribution_
virtual void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, const size_t inputBeginGlobalIndex, const uint32_t tofBin, const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield)
virtual void init(const edm::EventSetup &es)=0
Phase2TrackerDigitizerAlgorithm(const edm::ParameterSet &conf_common, const edm::ParameterSet &conf_specific, edm::ConsumesCollector iC)
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser_
virtual void induce_signal(const PSimHit &hit, const size_t hitIndex, const uint32_t tofBin, const Phase2TrackerGeomDetUnit *pixdet, const std::vector< DigitizerUtility::SignalPoint > &collection_points)
virtual void pixel_inefficiency_db(uint32_t detID)
std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Endcap_
virtual void digitize(const Phase2TrackerGeomDetUnit *pixdet, std::map< int, DigitizerUtility::DigiSimInfo > &digi_map, const TrackerTopology *tTopo)
std::map< int, DigitizerUtility::Amplitude, std::less< int > > signal_map_type
std::vector< edm::ParameterSet > Parameters
constexpr double c_inv
Definition: DetId.h:17
const SiPhase2OuterTrackerLorentzAngle * siPhase2OTLorentzAngle_
std::unique_ptr< CLHEP::RandGaussQ > smearedThreshold_Barrel_
virtual std::vector< DigitizerUtility::EnergyDepositUnit > primary_ionization(const PSimHit &hit) const
virtual void pixel_inefficiency(const SubdetEfficiencies &eff, const Phase2TrackerGeomDetUnit *pixdet, const TrackerTopology *tTopo)
virtual std::vector< DigitizerUtility::SignalPoint > drift(const PSimHit &hit, const Phase2TrackerGeomDetUnit *pixdet, const GlobalVector &bfield, const std::vector< DigitizerUtility::EnergyDepositUnit > &ionization_points) const
std::map< uint32_t, signal_map_type > signalMaps
virtual std::vector< float > fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments) const
int convertSignalToAdc(uint32_t detID, float signal_in_elec, float threshold)