CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiStripDigitizerAlgorithm.h
Go to the documentation of this file.
1 #ifndef SiStripDigitizerAlgorithm_h
2 #define SiStripDigitizerAlgorithm_h
3 
10 #include <memory>
11 #include <string>
12 
14 
30 #include "SiHitDigitizer.h"
41 
42 #include "TH1F.h"
43 
44 #include <iostream>
45 #include <fstream>
46 
47 class TrackerTopology;
48 
50 class StripDigiSimLink;
51 
52 namespace CLHEP {
53  class HepRandomEngine;
54 }
55 
57 public:
61  typedef std::map<int, float, std::less<int>> hit_map_type;
62  typedef float Amplitude;
63 
64  // Constructor
66 
67  // Destructor
69 
70  void initializeDetUnit(StripGeomDetUnit const* det, const edm::EventSetup& iSetup);
71 
72  void initializeEvent(const edm::EventSetup& iSetup);
73 
74  //run the algorithm to digitize a single det
75  void accumulateSimHits(const std::vector<PSimHit>::const_iterator inputBegin,
76  const std::vector<PSimHit>::const_iterator inputEnd,
77  size_t inputBeginGlobalIndex,
78  unsigned int tofBin,
79  const StripGeomDetUnit* stripdet,
80  const GlobalVector& bfield,
81  const TrackerTopology* tTopo,
82  CLHEP::HepRandomEngine*);
83 
84  void digitize(edm::DetSet<SiStripDigi>& outDigis,
85  edm::DetSet<SiStripRawDigi>& outRawDigis,
86  edm::DetSet<SiStripRawDigi>& outStripAmplitudes,
87  edm::DetSet<SiStripRawDigi>& outStripAmplitudesPostAPV,
88  edm::DetSet<SiStripRawDigi>& outStripAPVBaselines,
90  const StripGeomDetUnit* stripdet,
91  const SiStripGain&,
92  const SiStripThreshold&,
93  const SiStripNoises&,
94  const SiStripPedestals&,
95  bool simulateAPVInThisEvent,
97  std::vector<std::pair<int, std::bitset<6>>>& theAffectedAPVvector,
98  CLHEP::HepRandomEngine*,
99  const TrackerTopology* tTopo);
100 
102 
103  // ParticleDataTable
105  theSiHitDigitizer->setParticleDataTable(pardt);
106  pdt = pardt;
107  }
108 
109 private:
110  const double theThreshold;
111  const double cmnRMStib;
112  const double cmnRMStob;
113  const double cmnRMStid;
114  const double cmnRMStec;
116  const bool
117  makeDigiSimLinks_; //< Whether or not to create the association to sim truth collection. Set in configuration.
118  const bool peakMode;
119  const bool noise;
120  const bool RealPedestals;
121  const bool SingleStripNoise;
122  const bool CommonModeNoise;
123  const bool BaselineShift;
125 
126  const int theFedAlgo;
127  const bool zeroSuppression;
128  const double theElectronPerADC;
129 
130  const double theTOFCutForPeak;
132  const double tofCut;
133  const double cosmicShift;
134  const double inefficiency;
135  const double pedOffset;
136  const bool PreMixing_;
140 
143 
147 
148  const std::unique_ptr<SiHitDigitizer> theSiHitDigitizer;
149  const std::unique_ptr<SiPileUpSignals> theSiPileUpSignals;
150  const std::unique_ptr<const SiGaussianTailNoiseAdder> theSiNoiseAdder;
151  const std::unique_ptr<SiTrivialDigitalConverter> theSiDigitalConverter;
152  const std::unique_ptr<SiStripFedZeroSuppression> theSiZeroSuppress;
153 
154  // bad channels for each detector ID
155  std::map<unsigned int, std::vector<bool>> allBadChannels;
156  std::map<unsigned int, std::vector<bool>> allHIPChannels;
157  // first and last channel wit signal for each detector ID
158  std::map<unsigned int, size_t> firstChannelsWithSignal;
159  std::map<unsigned int, size_t> lastChannelsWithSignal;
160 
161  // ESHandles
163 
167  unsigned int trackID;
171  unsigned int tofBin; // Needed along with subDet to determine which PSimHit collection simHitGlobalIndex indexes
172  };
173 
174  typedef std::map<int, std::vector<AssociationInfo>> AssociationInfoForChannel;
175  typedef std::map<uint32_t, AssociationInfoForChannel> AssociationInfoForDetId;
178 
180 
181  std::ifstream APVProbaFile;
182  std::map<int, float> mapOfAPVprobabilities;
183  std::map<int, std::bitset<6>> SiStripTrackerAffectedAPVMap;
185 
187  const double apv_maxResponse_;
188  const double apv_rate_;
189  const double apv_mVPerQ_;
190  const double apv_fCPerElectron_;
191  unsigned int nTruePU_;
192 };
193 
194 #endif
const edm::ESGetToken< HepPDT::ParticleDataTable, PDTRecord > pdtToken_
const std::unique_ptr< SiPileUpSignals > theSiPileUpSignals
std::vector< SiStripDigi > DigitalVecType
edm::ESHandle< SiStripLorentzAngle > lorentzAngleHandle
std::map< unsigned int, std::vector< bool > > allHIPChannels
std::map< int, float > mapOfAPVprobabilities
HepPDT::ParticleDataTable ParticleDataTable
AssociationInfoForDetId associationInfoForDetId_
Structure that holds the information on the SimTrack contributions. Only filled if makeDigiSimLinks_ ...
SiDigitalConverter::DigitalRawVecType DigitalRawVecType
SiDigitalConverter::DigitalVecType DigitalVecType
void initializeEvent(const edm::EventSetup &iSetup)
void initializeDetUnit(StripGeomDetUnit const *det, const edm::EventSetup &iSetup)
const edm::ESGetToken< SiStripBadStrip, SiStripBadChannelRcd > deadChannelToken_
const std::unique_ptr< const SiGaussianTailNoiseAdder > theSiNoiseAdder
const std::unique_ptr< SiHitDigitizer > theSiHitDigitizer
std::map< unsigned int, size_t > firstChannelsWithSignal
std::map< int, Amplitude > SignalMapType
std::map< int, float, std::less< int > > hit_map_type
const std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
SiPileUpSignals::SignalMapType SignalMapType
const std::unique_ptr< SiStripFedZeroSuppression > theSiZeroSuppress
std::vector< SiStripRawDigi > DigitalRawVecType
void calculateInstlumiScale(PileupMixingContent *puInfo)
std::map< unsigned int, std::vector< bool > > allBadChannels
HepPDT::ParticleData ParticleData
void setParticleDataTable(const ParticleDataTable *pardt)
const edm::ESGetToken< SiStripLorentzAngle, SiStripLorentzAngleSimRcd > lorentzAngleToken_
std::map< int, std::vector< AssociationInfo > > AssociationInfoForChannel
SiStripDigitizerAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
std::map< unsigned int, size_t > lastChannelsWithSignal
void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, size_t inputBeginGlobalIndex, unsigned int tofBin, const StripGeomDetUnit *stripdet, const GlobalVector &bfield, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
std::map< uint32_t, AssociationInfoForChannel > AssociationInfoForDetId
std::map< int, std::bitset< 6 > > SiStripTrackerAffectedAPVMap
const ParticleDataTable * pdt
void digitize(edm::DetSet< SiStripDigi > &outDigis, edm::DetSet< SiStripRawDigi > &outRawDigis, edm::DetSet< SiStripRawDigi > &outStripAmplitudes, edm::DetSet< SiStripRawDigi > &outStripAmplitudesPostAPV, edm::DetSet< SiStripRawDigi > &outStripAPVBaselines, edm::DetSet< StripDigiSimLink > &outLink, const StripGeomDetUnit *stripdet, const SiStripGain &, const SiStripThreshold &, const SiStripNoises &, const SiStripPedestals &, bool simulateAPVInThisEvent, const SiStripApvSimulationParameters *, std::vector< std::pair< int, std::bitset< 6 >>> &theAffectedAPVvector, CLHEP::HepRandomEngine *, const TrackerTopology *tTopo)
size_t simHitGlobalIndex
The array index of the sim hit, but in the array for all crossings.