CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripDigitizerAlgorithm.cc
Go to the documentation of this file.
1 // File: SiStripDigitizerAlgorithm.cc
2 // Description: Steering class for digitization.
3 // Modified 15/May/2013 mark.grimes@bristol.ac.uk - Modified so that the digi-sim link has the correct
4 // index for the sim hits stored. It was previously always set to zero (I won't mention that it was
5 // me who originally wrote that).
6 // Modified on Feb 11, 2015: prolay.kumar.mal@cern.ch & Jean-Laurent.Agram@cern.ch
7 // Added/Modified the individual strip noise in zero suppression
8 // mode from the conditions DB; previously, the digitizer used to
9 // consider the noise value for individual strips inside a module from
10 // the central strip noise value.
12 #include <vector>
13 #include <algorithm>
14 #include <iostream>
31 #include "CLHEP/Random/RandFlat.h"
32 
34  lorentzAngleName(conf.getParameter<std::string>("LorentzAngle")),
35  theThreshold(conf.getParameter<double>("NoiseSigmaThreshold")),
36  cmnRMStib(conf.getParameter<double>("cmnRMStib")),
37  cmnRMStob(conf.getParameter<double>("cmnRMStob")),
38  cmnRMStid(conf.getParameter<double>("cmnRMStid")),
39  cmnRMStec(conf.getParameter<double>("cmnRMStec")),
40  APVSaturationProb(conf.getParameter<double>("APVSaturationProb")),
41  makeDigiSimLinks_(conf.getUntrackedParameter<bool>("makeDigiSimLinks", false)),
42  peakMode(conf.getParameter<bool>("APVpeakmode")),
43  noise(conf.getParameter<bool>("Noise")),
44  RealPedestals(conf.getParameter<bool>("RealPedestals")),
45  SingleStripNoise(conf.getParameter<bool>("SingleStripNoise")),
46  CommonModeNoise(conf.getParameter<bool>("CommonModeNoise")),
47  BaselineShift(conf.getParameter<bool>("BaselineShift")),
48  APVSaturationFromHIP(conf.getParameter<bool>("APVSaturationFromHIP")),
49  theFedAlgo(conf.getParameter<int>("FedAlgorithm")),
50  zeroSuppression(conf.getParameter<bool>("ZeroSuppression")),
51  theElectronPerADC(conf.getParameter<double>( peakMode ? "electronPerAdcPeak" : "electronPerAdcDec" )),
52  theTOFCutForPeak(conf.getParameter<double>("TOFCutForPeak")),
53  theTOFCutForDeconvolution(conf.getParameter<double>("TOFCutForDeconvolution")),
54  tofCut(peakMode ? theTOFCutForPeak : theTOFCutForDeconvolution),
55  cosmicShift(conf.getUntrackedParameter<double>("CosmicDelayShift")),
56  inefficiency(conf.getParameter<double>("Inefficiency")),
57  pedOffset((unsigned int)conf.getParameter<double>("PedestalsOffset")),
58  PreMixing_(conf.getParameter<bool>("PreMixingMode")),
59  theSiHitDigitizer(new SiHitDigitizer(conf)),
60  theSiPileUpSignals(new SiPileUpSignals()),
61  theSiNoiseAdder(new SiGaussianTailNoiseAdder(theThreshold)),
62  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC, PreMixing_)),
63  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo)) {
64 
65  if (peakMode) {
66  LogDebug("StripDigiInfo")<<"APVs running in peak mode (poor time resolution)";
67  } else {
68  LogDebug("StripDigiInfo")<<"APVs running in deconvolution mode (good time resolution)";
69  };
70  if(SingleStripNoise) LogDebug("SiStripDigitizerAlgorithm")<<" SingleStripNoise: ON";
71  else LogDebug("SiStripDigitizerAlgorithm")<<" SingleStripNoise: OFF";
72  if(CommonModeNoise) LogDebug("SiStripDigitizerAlgorithm")<<" CommonModeNoise: ON";
73  else LogDebug("SiStripDigitizerAlgorithm")<<" CommonModeNoise: OFF";
74 }
75 
77 }
78 
79 void
81  edm::ESHandle<SiStripBadStrip> deadChannelHandle;
82  iSetup.get<SiStripBadChannelRcd>().get(deadChannelHandle);
83 
84  unsigned int detId = det->geographicalId().rawId();
85  int numStrips = (det->specificTopology()).nstrips();
86 
87  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detId);
88  //storing the bad strip of the the module. the module is not removed but just signal put to 0
89  std::vector<bool>& badChannels = allBadChannels[detId];
90  badChannels.clear();
91  badChannels.insert(badChannels.begin(), numStrips, false);
92  for(SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
93  SiStripBadStrip::data fs = deadChannelHandle->decode(*it);
94  for(int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip) badChannels[strip] = true;
95  }
96  firstChannelsWithSignal[detId] = numStrips;
97  lastChannelsWithSignal[detId]= 0;
98 }
99 
100 void
102  theSiPileUpSignals->reset();
103  // This should be clear by after all calls to digitize(), but I might as well make sure
104  associationInfoForDetId_.clear();
105 
106  //get gain noise pedestal lorentzAngle from ES handle
108  iSetup.getData(pdt);
109  setParticleDataTable(&*pdt);
111 }
112 
113 // Run the algorithm for a given module
114 // ------------------------------------
115 
116 void
117 SiStripDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
118  std::vector<PSimHit>::const_iterator inputEnd,
119  size_t inputBeginGlobalIndex,
120  unsigned int tofBin,
121  const StripGeomDetUnit* det,
122  const GlobalVector& bfield,
123  const TrackerTopology *tTopo,
124  CLHEP::HepRandomEngine* engine) {
125  // produce SignalPoints for all SimHits in detector
126  unsigned int detID = det->geographicalId().rawId();
127  int numStrips = (det->specificTopology()).nstrips();
128 
129  std::vector<bool>& badChannels = allBadChannels[detID];
130  size_t thisFirstChannelWithSignal = numStrips;
131  size_t thisLastChannelWithSignal = 0;
132 
133  float langle = (lorentzAngleHandle.isValid()) ? lorentzAngleHandle->getLorentzAngle(detID) : 0.;
134 
135  std::vector<float> locAmpl(numStrips, 0.);
136 
137  // Loop over hits
138 
139  uint32_t detId = det->geographicalId().rawId();
140  // First: loop on the SimHits
141  if(CLHEP::RandFlat::shoot(engine) > inefficiency) {
142  AssociationInfoForChannel* pDetIDAssociationInfo; // I only need this if makeDigiSimLinks_ is true...
143  if( makeDigiSimLinks_ ) pDetIDAssociationInfo=&(associationInfoForDetId_[detId]); // ...so only search the map if that is the case
144  std::vector<float> previousLocalAmplitude; // Only used if makeDigiSimLinks_ is true. Needed to work out the change in amplitude.
145 
146  size_t simHitGlobalIndex=inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
147  for (std::vector<PSimHit>::const_iterator simHitIter = inputBegin; simHitIter != inputEnd; ++simHitIter, ++simHitGlobalIndex ) {
148  // skip hits not in this detector.
149  if((*simHitIter).detUnitId() != detId) {
150  continue;
151  }
152  // check TOF
153  if (std::fabs(simHitIter->tof() - cosmicShift - det->surface().toGlobal(simHitIter->localPosition()).mag()/30.) < tofCut && simHitIter->energyLoss()>0) {
154  if( makeDigiSimLinks_ ) previousLocalAmplitude=locAmpl; // Not needed except to make the sim link association.
155  size_t localFirstChannel = numStrips;
156  size_t localLastChannel = 0;
157  // process the hit
158  theSiHitDigitizer->processHit(&*simHitIter, *det, bfield, langle, locAmpl, localFirstChannel, localLastChannel, tTopo, engine);
159 
160  //APV Killer to simulate HIP effect
161  //------------------------------------------------------
162 
164  int pdg_id = simHitIter->particleType();
165  particle = pdt->particle(pdg_id);
166  if(particle != NULL){
167  float charge = particle->charge();
168  bool isHadron = particle->isHadron();
169  if(charge!=0 && isHadron){
170  if(CLHEP::RandFlat::shoot(engine) < APVSaturationProb){
171  int FirstAPV = localFirstChannel/128;
172  int LastAPV = localLastChannel/128;
173  //std::cout << "-------------------HIP--------------" << std::endl;
174  //std::cout << "Killing APVs " << FirstAPV << " - " <<LastAPV << " " << detID <<std::endl;
175  for(int strip = FirstAPV*128; strip < LastAPV*128 +128; ++strip) {
176  badChannels[strip] = true;
177  }
178  //doing like that I remove the signal information only after the
179  //stip that got the HIP but it remains the signal of the previous
180  //one. I'll make a further loop to remove all signal
181  }
182  }
183  }
184  }
185 
186 
187  if(thisFirstChannelWithSignal > localFirstChannel) thisFirstChannelWithSignal = localFirstChannel;
188  if(thisLastChannelWithSignal < localLastChannel) thisLastChannelWithSignal = localLastChannel;
189 
190  if( makeDigiSimLinks_ ) { // No need to do any of this if truth association was turned off in the configuration
191  for( size_t stripIndex=0; stripIndex<locAmpl.size(); ++stripIndex ) {
192  // Work out the amplitude from this SimHit from the difference of what it was before and what it is now
193  float signalFromThisSimHit=locAmpl[stripIndex]-previousLocalAmplitude[stripIndex];
194  if( signalFromThisSimHit!=0 ) { // If this SimHit had any contribution I need to record it.
195  auto& associationVector=(*pDetIDAssociationInfo)[stripIndex];
196  bool addNewEntry=true;
197  // Make sure the hit isn't in already. I've seen this a few times, it always seems to happen in pairs so I think
198  // it's something to do with the stereo strips.
199  for( auto& associationInfo : associationVector ) {
200  if( associationInfo.trackID==simHitIter->trackId() && associationInfo.eventID==simHitIter->eventId() ) {
201  // The hit is already in, so add this second contribution and move on
202  associationInfo.contributionToADC+=signalFromThisSimHit;
203  addNewEntry=false;
204  break;
205  }
206  } // end of loop over associationVector
207  // If the hit wasn't already in create a new association info structure.
208  if( addNewEntry ) associationVector.push_back( AssociationInfo{ simHitIter->trackId(), simHitIter->eventId(), signalFromThisSimHit, simHitGlobalIndex, tofBin } );
209  } // end of "if( signalFromThisSimHit!=0 )"
210  } // end of loop over locAmpl strips
211  } // end of "if( makeDigiSimLinks_ )"
212  } // end of TOF check
213  } // end for
214  }
215  theSiPileUpSignals->add(detID, locAmpl, thisFirstChannelWithSignal, thisLastChannelWithSignal);
216 
217  if(firstChannelsWithSignal[detID] > thisFirstChannelWithSignal) firstChannelsWithSignal[detID] = thisFirstChannelWithSignal;
218  if(lastChannelsWithSignal[detID] < thisLastChannelWithSignal) lastChannelsWithSignal[detID] = thisLastChannelWithSignal;
219 }
220 
221 void
223  edm::DetSet<SiStripDigi>& outdigi,
224  edm::DetSet<SiStripRawDigi>& outrawdigi,
226  const StripGeomDetUnit *det,
227  edm::ESHandle<SiStripGain> & gainHandle,
228  edm::ESHandle<SiStripThreshold> & thresholdHandle,
229  edm::ESHandle<SiStripNoises> & noiseHandle,
230  edm::ESHandle<SiStripPedestals> & pedestalHandle,
231  CLHEP::HepRandomEngine* engine) {
232  unsigned int detID = det->geographicalId().rawId();
233  int numStrips = (det->specificTopology()).nstrips();
234 
235  const SiPileUpSignals::SignalMapType* theSignal(theSiPileUpSignals->getSignal(detID));
236 
237  std::vector<float> detAmpl(numStrips, 0.);
238  if(theSignal) {
239  for(const auto& amp : *theSignal) {
240  detAmpl[amp.first] = amp.second;
241  }
242  }
243 
244  //removing signal from the dead (and HIP effected) strips
245  std::vector<bool>& badChannels = allBadChannels[detID];
246  for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] = 0.;
247 
248  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
249  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
250  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
251 
252 // -----------------------------------------------------------
253 
254  auto& firstChannelWithSignal = firstChannelsWithSignal[detID];
255  auto& lastChannelWithSignal = lastChannelsWithSignal[detID];
256  auto iAssociationInfoByChannel=associationInfoForDetId_.find(detID); // Use an iterator so that I can easily remove it once finished
257 
258  if(zeroSuppression){
259 
260  //Adding the strip noise
261  //------------------------------------------------------
262  if(noise){
263 
264  if(SingleStripNoise){
265  std::vector<float> noiseRMSv;
266  noiseRMSv.clear();
267  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
268  for(int strip=0; strip< numStrips; ++strip){
269  if(!badChannels[strip]){
270  float gainValue = gainHandle->getStripGain(strip, detGainRange);
271  noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC/gainValue;
272  //std::cout<<"<SiStripDigitizerAlgorithm::digitize>: gainValue: "<<gainValue<<"\tnoiseRMSv["<<strip<<"]: "<<noiseRMSv[strip]<<std::endl;
273  }
274  }
275  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
276  } else {
277  int RefStrip = int(numStrips/2.);
278  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
279  RefStrip++;
280  }
281  if(RefStrip<numStrips){
282  float RefgainValue = gainHandle->getStripGain(RefStrip, detGainRange);
283  float RefnoiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC/RefgainValue;
284 
285  theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,RefnoiseRMS, engine);
286  //std::cout<<"<SiStripDigitizerAlgorithm::digitize>: RefgainValue: "<<RefgainValue<<"\tRefnoiseRMS: "<<RefnoiseRMS<<std::endl;
287  }
288  }
289  }//if noise
290 
291  DigitalVecType digis;
292  theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
293  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
294  // will be empty and the iterator will always equal associationInfoForDetId_.end().
295  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
296  for( const auto& iDigi : digis ) {
297  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
298  const std::vector<AssociationInfo>& associationInfo=associationInfoByChannel[iDigi.channel()];
299 
300  // Need to find the total from all sim hits, because this might not be the same as the total
301  // digitised due to noise or whatever.
302  float totalSimADC=0;
303  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
304  // Now I know that I can loop again and create the links
305  for( const auto& iAssociationInfo : associationInfo ) {
306  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
307  // expected counting to start at 1, not 0. Now changed.
308  outLink.push_back( StripDigiSimLink( iDigi.channel(), iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex, iAssociationInfo.tofBin, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
309  } // end of loop over associationInfo
310  } // end of loop over the digis
311  } // end of check that iAssociationInfoByChannel is a valid iterator
312  outdigi.data = digis;
313  }//if zeroSuppression
314 
315  if(!zeroSuppression){
316  //if(noise){
317  // the constant pedestal offset is needed because
318  // negative adc counts are not allowed in case
319  // Pedestal and CMN subtraction is performed.
320  // The pedestal value read from the conditions
321  // is pedValue and after the pedestal subtraction
322  // the baseline is zero. The Common Mode Noise
323  // is not subtracted from the negative adc counts
324  // channels. Adding pedOffset the baseline is set
325  // to pedOffset after pedestal subtraction and CMN
326  // is subtracted to all the channels since none of
327  // them has negative adc value. The pedOffset is
328  // treated as a constant component in the CMN
329  // estimation and subtracted as CMN.
330 
331 
332  //calculating the charge deposited on each APV and subtracting the shift
333  //------------------------------------------------------
334  if(BaselineShift){
335  theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
336  }
337 
338  //Adding the strip noise
339  //------------------------------------------------------
340  if(noise){
341  std::vector<float> noiseRMSv;
342  noiseRMSv.clear();
343  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
344 
345  if(SingleStripNoise){
346  for(int strip=0; strip< numStrips; ++strip){
347  if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
348  }
349 
350  } else {
351  int RefStrip = 0; //int(numStrips/2.);
352  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
353  RefStrip++;
354  }
355  if(RefStrip<numStrips){
356  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
357  for(int strip=0; strip< numStrips; ++strip){
358  if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
359  }
360  }
361  }
362 
363  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
364  }
365 
366  //adding the CMN
367  //------------------------------------------------------
368  if(CommonModeNoise){
369  float cmnRMS = 0.;
370  DetId detId(detID);
371  uint32_t SubDet = detId.subdetId();
372  if(SubDet==3){
373  cmnRMS = cmnRMStib;
374  }else if(SubDet==4){
375  cmnRMS = cmnRMStid;
376  }else if(SubDet==5){
377  cmnRMS = cmnRMStob;
378  }else if(SubDet==6){
379  cmnRMS = cmnRMStec;
380  }
381  cmnRMS *= theElectronPerADC;
382  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels, engine);
383  }
384 
385 
386  //Adding the pedestals
387  //------------------------------------------------------
388 
389  std::vector<float> vPeds;
390  vPeds.clear();
391  vPeds.insert(vPeds.begin(),numStrips,0.);
392 
393  if(RealPedestals){
394  for(int strip=0; strip< numStrips; ++strip){
395  if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
396  }
397  } else {
398  for(int strip=0; strip< numStrips; ++strip){
399  if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
400  }
401  }
402 
403  theSiNoiseAdder->addPedestals(detAmpl, vPeds);
404 
405 
406  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
407  // edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
408  //}else{
409 
410  DigitalRawVecType rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
411 
412  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
413  // will be empty and the iterator will always equal associationInfoForDetId_.end().
414  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
415  // N.B. For the raw digis the channel is inferred from the position in the vector.
416  // I'VE NOT TESTED THIS YET!!!!!
417  // ToDo Test this properly.
418  for( size_t channel=0; channel<rawdigis.size(); ++channel ) {
419  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
420  const auto iAssociationInfo=associationInfoByChannel.find(channel);
421  if( iAssociationInfo==associationInfoByChannel.end() ) continue; // Skip if there is no sim information for this channel (i.e. it's all noise)
422  const std::vector<AssociationInfo>& associationInfo=iAssociationInfo->second;
423 
424  // Need to find the total from all sim hits, because this might not be the same as the total
425  // digitised due to noise or whatever.
426  float totalSimADC=0;
427  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
428  // Now I know that I can loop again and create the links
429  for( const auto& iAssociationInfo : associationInfo ) {
430  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
431  // expected counting to start at 1, not 0. Now changed.
432  outLink.push_back( StripDigiSimLink( channel, iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex, iAssociationInfo.tofBin, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
433  } // end of loop over associationInfo
434  } // end of loop over the digis
435  } // end of check that iAssociationInfoByChannel is a valid iterator
436 
437  outrawdigi.data = rawdigis;
438 
439  //}
440  }
441 
442  // Now that I've finished with this entry in the map of associations, I can remove it.
443  // Note that there might not be an association if the ADC reading is from noise in which
444  // case associationIsValid will be false.
445  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) associationInfoForDetId_.erase(iAssociationInfoByChannel);
446 }
#define LogDebug(id)
unsigned short range
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
const std::unique_ptr< SiPileUpSignals > theSiPileUpSignals
edm::ESHandle< SiStripLorentzAngle > lorentzAngleHandle
void push_back(const T &t)
Definition: DetSet.h:68
AssociationInfoForDetId associationInfoForDetId_
Structure that holds the information on the SimTrack contributions. Only filled if makeDigiSimLinks_ ...
SiDigitalConverter::DigitalRawVecType DigitalRawVecType
SiDigitalConverter::DigitalVecType DigitalVecType
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< unsigned int >::const_iterator ContainerIterator
void initializeEvent(const edm::EventSetup &iSetup)
#define NULL
Definition: scimark2.h:8
std::pair< ContainerIterator, ContainerIterator > Range
void initializeDetUnit(StripGeomDetUnit const *det, const edm::EventSetup &iSetup)
void digitize(edm::DetSet< SiStripDigi > &outDigis, edm::DetSet< SiStripRawDigi > &outRawDigis, edm::DetSet< StripDigiSimLink > &outLink, const StripGeomDetUnit *stripdet, edm::ESHandle< SiStripGain > &, edm::ESHandle< SiStripThreshold > &, edm::ESHandle< SiStripNoises > &, edm::ESHandle< SiStripPedestals > &, CLHEP::HepRandomEngine *)
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
const std::unique_ptr< const SiGaussianTailNoiseAdder > theSiNoiseAdder
const std::unique_ptr< SiHitDigitizer > theSiHitDigitizer
void getData(T &iHolder) const
Definition: EventSetup.h:79
std::map< unsigned int, size_t > firstChannelsWithSignal
std::map< int, Amplitude > SignalMapType
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
const std::unique_ptr< SiStripFedZeroSuppression > theSiZeroSuppress
std::pair< ContainerIterator, ContainerIterator > Range
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:77
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::map< unsigned int, std::vector< bool > > allBadChannels
void setParticleDataTable(const ParticleDataTable *pardt)
Definition: DetId.h:18
unsigned short firstStrip
std::map< int, std::vector< AssociationInfo > > AssociationInfoForChannel
const T & get() const
Definition: EventSetup.h:56
SiStripDigitizerAlgorithm(const edm::ParameterSet &conf)
collection_type data
Definition: DetSet.h:78
std::pair< ContainerIterator, ContainerIterator > Range
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 *)
volatile std::atomic< bool > shutdown_flag false
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:48
bool isValid() const
Definition: ESHandle.h:47
const ParticleDataTable * pdt