CMS 3D CMS Logo

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 
33 #include <cstring>
34 #include <sstream>
35 
36 #include <boost/algorithm/string.hpp>
37 
39  : lorentzAngleName(conf.getParameter<std::string>("LorentzAngle")),
40  theThreshold(conf.getParameter<double>("NoiseSigmaThreshold")),
41  cmnRMStib(conf.getParameter<double>("cmnRMStib")),
42  cmnRMStob(conf.getParameter<double>("cmnRMStob")),
43  cmnRMStid(conf.getParameter<double>("cmnRMStid")),
44  cmnRMStec(conf.getParameter<double>("cmnRMStec")),
45  APVSaturationProbScaling_(conf.getParameter<double>("APVSaturationProbScaling")),
46  makeDigiSimLinks_(conf.getUntrackedParameter<bool>("makeDigiSimLinks", false)),
47  peakMode(conf.getParameter<bool>("APVpeakmode")),
48  noise(conf.getParameter<bool>("Noise")),
49  RealPedestals(conf.getParameter<bool>("RealPedestals")),
50  SingleStripNoise(conf.getParameter<bool>("SingleStripNoise")),
51  CommonModeNoise(conf.getParameter<bool>("CommonModeNoise")),
52  BaselineShift(conf.getParameter<bool>("BaselineShift")),
53  APVSaturationFromHIP(conf.getParameter<bool>("APVSaturationFromHIP")),
54  theFedAlgo(conf.getParameter<int>("FedAlgorithm")),
55  zeroSuppression(conf.getParameter<bool>("ZeroSuppression")),
56  theElectronPerADC(conf.getParameter<double>(peakMode ? "electronPerAdcPeak" : "electronPerAdcDec")),
57  theTOFCutForPeak(conf.getParameter<double>("TOFCutForPeak")),
58  theTOFCutForDeconvolution(conf.getParameter<double>("TOFCutForDeconvolution")),
59  tofCut(peakMode ? theTOFCutForPeak : theTOFCutForDeconvolution),
60  cosmicShift(conf.getUntrackedParameter<double>("CosmicDelayShift")),
61  inefficiency(conf.getParameter<double>("Inefficiency")),
62  pedOffset((unsigned int)conf.getParameter<double>("PedestalsOffset")),
63  PreMixing_(conf.getParameter<bool>("PreMixingMode")),
64  theSiHitDigitizer(new SiHitDigitizer(conf)),
65  theSiPileUpSignals(new SiPileUpSignals()),
66  theSiNoiseAdder(new SiGaussianTailNoiseAdder(theThreshold)),
67  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC, PreMixing_)),
68  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo)),
69  APVProbabilityFile(conf.getParameter<edm::FileInPath>("APVProbabilityFile")) {
70  if (peakMode) {
71  LogDebug("StripDigiInfo") << "APVs running in peak mode (poor time resolution)";
72  } else {
73  LogDebug("StripDigiInfo") << "APVs running in deconvolution mode (good time resolution)";
74  };
75  if (SingleStripNoise)
76  LogDebug("SiStripDigitizerAlgorithm") << " SingleStripNoise: ON";
77  else
78  LogDebug("SiStripDigitizerAlgorithm") << " SingleStripNoise: OFF";
79  if (CommonModeNoise)
80  LogDebug("SiStripDigitizerAlgorithm") << " CommonModeNoise: ON";
81  else
82  LogDebug("SiStripDigitizerAlgorithm") << " CommonModeNoise: OFF";
84  throw cms::Exception("PreMixing does not work with HIP loss simulation yet");
87  APVProbaFile.open((APVProbabilityFile.fullPath()).c_str());
88  if (APVProbaFile.is_open()) {
89  while (getline(APVProbaFile, line)) {
90  std::vector<std::string> strs;
91  boost::split(strs, line, boost::is_any_of(" "));
92  if (strs.size() == 2) {
93  mapOfAPVprobabilities[std::stoi(strs.at(0))] = std::stof(strs.at(1));
94  }
95  }
96  APVProbaFile.close();
97  } else
98  throw cms::Exception("MissingInput") << "It seems that the APV probability list is missing\n";
99  }
100 }
101 
103 
105  edm::ESHandle<SiStripBadStrip> deadChannelHandle;
106  iSetup.get<SiStripBadChannelRcd>().get(deadChannelHandle);
107 
108  unsigned int detId = det->geographicalId().rawId();
109  int numStrips = (det->specificTopology()).nstrips();
110 
111  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detId);
112  //storing the bad strip of the the module. the module is not removed but just signal put to 0
113  std::vector<bool>& badChannels = allBadChannels[detId];
114  badChannels.clear();
115  badChannels.insert(badChannels.begin(), numStrips, false);
116  for (SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
117  SiStripBadStrip::data fs = deadChannelHandle->decode(*it);
118  for (int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip) {
119  badChannels[strip] = true;
120  }
121  }
122  firstChannelsWithSignal[detId] = numStrips;
123  lastChannelsWithSignal[detId] = 0;
124 
125  // if(APVSaturationFromHIP){
126  // std::bitset<6> &bs=SiStripTrackerAffectedAPVMap[detId];
127  // if(bs.any())theAffectedAPVvector.push_back(std::make_pair(detId,bs));
128  //}
129 }
130 
132  theSiPileUpSignals->reset();
133  // This should be clear by after all calls to digitize(), but I might as well make sure
134  associationInfoForDetId_.clear();
135 
136  APVSaturationProb_ = APVSaturationProbScaling_; // reset probability
138  FirstLumiCalc_ = true;
139  FirstDigitize_ = true;
140 
141  //get gain noise pedestal lorentzAngle from ES handle
143  iSetup.getData(pdt);
144  setParticleDataTable(&*pdt);
146 }
147 
148 // Run the algorithm for a given module
149 // ------------------------------------
150 
151 void SiStripDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
152  std::vector<PSimHit>::const_iterator inputEnd,
153  size_t inputBeginGlobalIndex,
154  unsigned int tofBin,
155  const StripGeomDetUnit* det,
156  const GlobalVector& bfield,
157  const TrackerTopology* tTopo,
158  CLHEP::HepRandomEngine* engine) {
159  // produce SignalPoints for all SimHits in detector
160  unsigned int detID = det->geographicalId().rawId();
161  int numStrips = (det->specificTopology()).nstrips();
162 
163  size_t thisFirstChannelWithSignal = numStrips;
164  size_t thisLastChannelWithSignal = 0;
165 
166  float langle = (lorentzAngleHandle.isValid()) ? lorentzAngleHandle->getLorentzAngle(detID) : 0.;
167 
168  std::vector<float> locAmpl(numStrips, 0.);
169 
170  // Loop over hits
171 
172  uint32_t detId = det->geographicalId().rawId();
173  // First: loop on the SimHits
174  if (CLHEP::RandFlat::shoot(engine) > inefficiency) {
175  AssociationInfoForChannel* pDetIDAssociationInfo; // I only need this if makeDigiSimLinks_ is true...
176  if (makeDigiSimLinks_)
177  pDetIDAssociationInfo = &(associationInfoForDetId_[detId]); // ...so only search the map if that is the case
178  std::vector<float>
179  previousLocalAmplitude; // Only used if makeDigiSimLinks_ is true. Needed to work out the change in amplitude.
180 
181  size_t simHitGlobalIndex = inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
182  for (std::vector<PSimHit>::const_iterator simHitIter = inputBegin; simHitIter != inputEnd;
183  ++simHitIter, ++simHitGlobalIndex) {
184  // skip hits not in this detector.
185  if ((*simHitIter).detUnitId() != detId) {
186  continue;
187  }
188  // check TOF
189  if (std::fabs(simHitIter->tof() - cosmicShift -
190  det->surface().toGlobal(simHitIter->localPosition()).mag() / 30.) < tofCut &&
191  simHitIter->energyLoss() > 0) {
192  if (makeDigiSimLinks_)
193  previousLocalAmplitude = locAmpl; // Not needed except to make the sim link association.
194  size_t localFirstChannel = numStrips;
195  size_t localLastChannel = 0;
196  // process the hit
197  theSiHitDigitizer->processHit(
198  &*simHitIter, *det, bfield, langle, locAmpl, localFirstChannel, localLastChannel, tTopo, engine);
199 
200  if (thisFirstChannelWithSignal > localFirstChannel)
201  thisFirstChannelWithSignal = localFirstChannel;
202  if (thisLastChannelWithSignal < localLastChannel)
203  thisLastChannelWithSignal = localLastChannel;
204 
205  if (makeDigiSimLinks_) { // No need to do any of this if truth association was turned off in the configuration
206  for (size_t stripIndex = 0; stripIndex < locAmpl.size(); ++stripIndex) {
207  // Work out the amplitude from this SimHit from the difference of what it was before and what it is now
208  float signalFromThisSimHit = locAmpl[stripIndex] - previousLocalAmplitude[stripIndex];
209  if (signalFromThisSimHit != 0) { // If this SimHit had any contribution I need to record it.
210  auto& associationVector = (*pDetIDAssociationInfo)[stripIndex];
211  bool addNewEntry = true;
212  // 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
213  // it's something to do with the stereo strips.
214  for (auto& associationInfo : associationVector) {
215  if (associationInfo.trackID == simHitIter->trackId() &&
216  associationInfo.eventID == simHitIter->eventId()) {
217  // The hit is already in, so add this second contribution and move on
218  associationInfo.contributionToADC += signalFromThisSimHit;
219  addNewEntry = false;
220  break;
221  }
222  } // end of loop over associationVector
223  // If the hit wasn't already in create a new association info structure.
224  if (addNewEntry)
225  associationVector.push_back(AssociationInfo{
226  simHitIter->trackId(), simHitIter->eventId(), signalFromThisSimHit, simHitGlobalIndex, tofBin});
227  } // end of "if( signalFromThisSimHit!=0 )"
228  } // end of loop over locAmpl strips
229  } // end of "if( makeDigiSimLinks_ )"
230  } // end of TOF check
231  } // end for
232  }
233  theSiPileUpSignals->add(detID, locAmpl, thisFirstChannelWithSignal, thisLastChannelWithSignal);
234 
235  if (firstChannelsWithSignal[detID] > thisFirstChannelWithSignal)
236  firstChannelsWithSignal[detID] = thisFirstChannelWithSignal;
237  if (lastChannelsWithSignal[detID] < thisLastChannelWithSignal)
238  lastChannelsWithSignal[detID] = thisLastChannelWithSignal;
239 }
240 
241 //============================================================================
243  //Instlumi scalefactor calculating for dynamic inefficiency
244 
245  if (puInfo && FirstLumiCalc_) {
246  const std::vector<int>& bunchCrossing = puInfo->getMix_bunchCrossing();
247  const std::vector<float>& TrueInteractionList = puInfo->getMix_TrueInteractions();
248  const int bunchSpacing = puInfo->getMix_bunchSpacing();
249 
250  double RevFreq = 11245.;
251  double minBXsec = 70.0E-27; // use 70mb as an approximation
252  double Bunch = 2100.; // 2016 value
253  if (bunchSpacing == 50)
254  Bunch = Bunch / 2.;
255 
256  int pui = 0, p = 0;
257  std::vector<int>::const_iterator pu;
258  std::vector<int>::const_iterator pu0 = bunchCrossing.end();
259 
260  for (pu = bunchCrossing.begin(); pu != bunchCrossing.end(); ++pu) {
261  if (*pu == 0) {
262  pu0 = pu;
263  p = pui;
264  }
265  pui++;
266  }
267  if (pu0 != bunchCrossing.end()) { // found the in-time interaction
268  double Tintr = TrueInteractionList.at(p);
269  double instLumi = Bunch * Tintr * RevFreq / minBXsec;
270  APVSaturationProb_ = instLumi / 6.0E33;
271  }
272  FirstLumiCalc_ = false;
273  }
274 }
275 
276 //============================================================================
277 
279  edm::DetSet<SiStripRawDigi>& outrawdigi,
281  const StripGeomDetUnit* det,
282  edm::ESHandle<SiStripGain>& gainHandle,
283  edm::ESHandle<SiStripThreshold>& thresholdHandle,
284  edm::ESHandle<SiStripNoises>& noiseHandle,
285  edm::ESHandle<SiStripPedestals>& pedestalHandle,
286  std::vector<std::pair<int, std::bitset<6>>>& theAffectedAPVvector,
287  CLHEP::HepRandomEngine* engine) {
288  unsigned int detID = det->geographicalId().rawId();
289  int numStrips = (det->specificTopology()).nstrips();
290 
291  const SiPileUpSignals::SignalMapType* theSignal(theSiPileUpSignals->getSignal(detID));
292 
293  std::vector<float> detAmpl(numStrips, 0.);
294  if (theSignal) {
295  for (const auto& amp : *theSignal) {
296  detAmpl[amp.first] = amp.second;
297  }
298  }
299 
300  //removing signal from the dead (and HIP effected) strips
301  std::vector<bool>& badChannels = allBadChannels[detID];
302  for (int strip = 0; strip < numStrips; ++strip) {
303  if (badChannels[strip]) {
304  detAmpl[strip] = 0.;
305  }
306  }
307 
308  if (APVSaturationFromHIP) {
309  //Implementation of the proper charge scaling function. Need consider resaturation effect:
310  //The probability map gives the probability that at least one HIP happened during the last N bunch crossings (cfr APV recovery time).
311  //The impact on the charge depends on the clostest HIP occurance (in terms of bunch crossing).
312  //The function discribing the APV recovery is therefore the weighted average function which takes into account all possibilities of HIP occurances across the last bx's.
313 
314  // do this step here because we now have access to luminosity information
315  if (FirstDigitize_) {
316  for (std::map<int, float>::iterator iter = mapOfAPVprobabilities.begin(); iter != mapOfAPVprobabilities.end();
317  ++iter) {
318  std::bitset<6> bs;
319  for (int Napv = 0; Napv < 6; Napv++) {
320  float cursor = CLHEP::RandFlat::shoot(engine);
321  bs[Napv] = cursor < iter->second * APVSaturationProb_
322  ? true
323  : false; //APVSaturationProb has been scaled by PU luminosity
324  }
325  SiStripTrackerAffectedAPVMap[iter->first] = bs;
326  }
327 
329  bool HasAtleastOneAffectedAPV = false;
330  while (!HasAtleastOneAffectedAPV) {
331  for (int bx = floor(300.0 / 25.0); bx > 0; bx--) { //Reminder: make these numbers not hard coded!!
332  float temp = CLHEP::RandFlat::shoot(engine) < 0.5 ? 1 : 0;
333  if (temp == 1 && bx < NumberOfBxBetweenHIPandEvent) {
335  HasAtleastOneAffectedAPV = true;
336  }
337  }
338  }
339 
340  FirstDigitize_ = false;
341  }
342 
343  std::bitset<6>& bs = SiStripTrackerAffectedAPVMap[detID];
344 
345  if (bs.any()) {
346  // store this information so it can be saved to the event later
347  theAffectedAPVvector.push_back(std::make_pair(detID, bs));
348 
349  if (!PreMixing_) {
350  // Here below is the scaling function which describes the evolution of the baseline (i.e. how the charge is suppressed).
351  // This must be replaced as soon as we have a proper modeling of the baseline evolution from VR runs
352  float Shift =
353  1 - NumberOfBxBetweenHIPandEvent / floor(300.0 / 25.0); //Reminder: make these numbers not hardcoded!!
354  float randomX = CLHEP::RandFlat::shoot(engine);
355  float scalingValue = (randomX - Shift) * 10.0 / 7.0 - 3.0 / 7.0;
356 
357  for (int strip = 0; strip < numStrips; ++strip) {
358  if (!badChannels[strip] && bs[strip / 128] == 1) {
359  detAmpl[strip] *= scalingValue > 0 ? scalingValue : 0.0;
360  }
361  }
362  }
363  }
364  }
365 
366  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
367  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
368  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
369 
370  // -----------------------------------------------------------
371 
372  auto& firstChannelWithSignal = firstChannelsWithSignal[detID];
373  auto& lastChannelWithSignal = lastChannelsWithSignal[detID];
374  auto iAssociationInfoByChannel =
375  associationInfoForDetId_.find(detID); // Use an iterator so that I can easily remove it once finished
376 
377  if (zeroSuppression) {
378  //Adding the strip noise
379  //------------------------------------------------------
380  if (noise) {
381  if (SingleStripNoise) {
382  // std::cout<<"In SSN, detId="<<detID<<std::endl;
383  std::vector<float> noiseRMSv;
384  noiseRMSv.clear();
385  noiseRMSv.insert(noiseRMSv.begin(), numStrips, 0.);
386  for (int strip = 0; strip < numStrips; ++strip) {
387  if (!badChannels[strip]) {
388  float gainValue = gainHandle->getStripGain(strip, detGainRange);
389  noiseRMSv[strip] = (noiseHandle->getNoise(strip, detNoiseRange)) * theElectronPerADC / gainValue;
390  //std::cout<<"<SiStripDigitizerAlgorithm::digitize>: gainValue: "<<gainValue<<"\tnoiseRMSv["<<strip<<"]: "<<noiseRMSv[strip]<<std::endl;
391  }
392  }
393  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
394  } else {
395  int RefStrip = int(numStrips / 2.);
396  while (RefStrip < numStrips &&
397  badChannels[RefStrip]) { //if the refstrip is bad, I move up to when I don't find it
398  RefStrip++;
399  }
400  if (RefStrip < numStrips) {
401  float RefgainValue = gainHandle->getStripGain(RefStrip, detGainRange);
402  float RefnoiseRMS = noiseHandle->getNoise(RefStrip, detNoiseRange) * theElectronPerADC / RefgainValue;
403 
404  theSiNoiseAdder->addNoise(
405  detAmpl, firstChannelWithSignal, lastChannelWithSignal, numStrips, RefnoiseRMS, engine);
406  //std::cout<<"<SiStripDigitizerAlgorithm::digitize>: RefgainValue: "<<RefgainValue<<"\tRefnoiseRMS: "<<RefnoiseRMS<<std::endl;
407  }
408  }
409  } //if noise
410 
411  DigitalVecType digis;
412  theSiZeroSuppress->suppress(
413  theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID, noiseHandle, thresholdHandle);
414  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
415  // will be empty and the iterator will always equal associationInfoForDetId_.end().
416  if (iAssociationInfoByChannel !=
417  associationInfoForDetId_.end()) { // make sure the readings for this DetID aren't completely from noise
418  for (const auto& iDigi : digis) {
419  auto& associationInfoByChannel = iAssociationInfoByChannel->second;
420  const std::vector<AssociationInfo>& associationInfo = associationInfoByChannel[iDigi.channel()];
421 
422  // Need to find the total from all sim hits, because this might not be the same as the total
423  // digitised due to noise or whatever.
424  float totalSimADC = 0;
425  for (const auto& iAssociationInfo : associationInfo)
426  totalSimADC += iAssociationInfo.contributionToADC;
427  // Now I know that I can loop again and create the links
428  for (const auto& iAssociationInfo : associationInfo) {
429  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
430  // expected counting to start at 1, not 0. Now changed.
431  outLink.push_back(StripDigiSimLink(iDigi.channel(),
432  iAssociationInfo.trackID,
433  iAssociationInfo.simHitGlobalIndex,
434  iAssociationInfo.tofBin,
435  iAssociationInfo.eventID,
436  iAssociationInfo.contributionToADC / totalSimADC));
437  } // end of loop over associationInfo
438  } // end of loop over the digis
439  } // end of check that iAssociationInfoByChannel is a valid iterator
440  outdigi.data = digis;
441  } //if zeroSuppression
442 
443  if (!zeroSuppression) {
444  //if(noise){
445  // the constant pedestal offset is needed because
446  // negative adc counts are not allowed in case
447  // Pedestal and CMN subtraction is performed.
448  // The pedestal value read from the conditions
449  // is pedValue and after the pedestal subtraction
450  // the baseline is zero. The Common Mode Noise
451  // is not subtracted from the negative adc counts
452  // channels. Adding pedOffset the baseline is set
453  // to pedOffset after pedestal subtraction and CMN
454  // is subtracted to all the channels since none of
455  // them has negative adc value. The pedOffset is
456  // treated as a constant component in the CMN
457  // estimation and subtracted as CMN.
458 
459  //calculating the charge deposited on each APV and subtracting the shift
460  //------------------------------------------------------
461  if (BaselineShift) {
462  theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
463  }
464 
465  //Adding the strip noise
466  //------------------------------------------------------
467  if (noise) {
468  std::vector<float> noiseRMSv;
469  noiseRMSv.clear();
470  noiseRMSv.insert(noiseRMSv.begin(), numStrips, 0.);
471 
472  if (SingleStripNoise) {
473  for (int strip = 0; strip < numStrips; ++strip) {
474  if (!badChannels[strip])
475  noiseRMSv[strip] = (noiseHandle->getNoise(strip, detNoiseRange)) * theElectronPerADC;
476  }
477 
478  } else {
479  int RefStrip = 0; //int(numStrips/2.);
480  while (RefStrip < numStrips &&
481  badChannels[RefStrip]) { //if the refstrip is bad, I move up to when I don't find it
482  RefStrip++;
483  }
484  if (RefStrip < numStrips) {
485  float noiseRMS = noiseHandle->getNoise(RefStrip, detNoiseRange) * theElectronPerADC;
486  for (int strip = 0; strip < numStrips; ++strip) {
487  if (!badChannels[strip])
488  noiseRMSv[strip] = noiseRMS;
489  }
490  }
491  }
492 
493  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
494  }
495 
496  //adding the CMN
497  //------------------------------------------------------
498  if (CommonModeNoise) {
499  float cmnRMS = 0.;
500  DetId detId(detID);
501  uint32_t SubDet = detId.subdetId();
502  if (SubDet == 3) {
503  cmnRMS = cmnRMStib;
504  } else if (SubDet == 4) {
505  cmnRMS = cmnRMStid;
506  } else if (SubDet == 5) {
507  cmnRMS = cmnRMStob;
508  } else if (SubDet == 6) {
509  cmnRMS = cmnRMStec;
510  }
511  cmnRMS *= theElectronPerADC;
512  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels, engine);
513  }
514 
515  //Adding the pedestals
516  //------------------------------------------------------
517 
518  std::vector<float> vPeds;
519  vPeds.clear();
520  vPeds.insert(vPeds.begin(), numStrips, 0.);
521 
522  if (RealPedestals) {
523  for (int strip = 0; strip < numStrips; ++strip) {
524  if (!badChannels[strip])
525  vPeds[strip] = (pedestalHandle->getPed(strip, detPedestalRange) + pedOffset) * theElectronPerADC;
526  }
527  } else {
528  for (int strip = 0; strip < numStrips; ++strip) {
529  if (!badChannels[strip])
530  vPeds[strip] = pedOffset * theElectronPerADC;
531  }
532  }
533 
534  theSiNoiseAdder->addPedestals(detAmpl, vPeds);
535 
536  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
537  // edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
538  //}else{
539 
540  DigitalRawVecType rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
541 
542  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
543  // will be empty and the iterator will always equal associationInfoForDetId_.end().
544  if (iAssociationInfoByChannel !=
545  associationInfoForDetId_.end()) { // make sure the readings for this DetID aren't completely from noise
546  // N.B. For the raw digis the channel is inferred from the position in the vector.
547  // I'VE NOT TESTED THIS YET!!!!!
548  // ToDo Test this properly.
549  for (size_t channel = 0; channel < rawdigis.size(); ++channel) {
550  auto& associationInfoByChannel = iAssociationInfoByChannel->second;
551  const auto iAssociationInfo = associationInfoByChannel.find(channel);
552  if (iAssociationInfo == associationInfoByChannel.end())
553  continue; // Skip if there is no sim information for this channel (i.e. it's all noise)
554  const std::vector<AssociationInfo>& associationInfo = iAssociationInfo->second;
555 
556  // Need to find the total from all sim hits, because this might not be the same as the total
557  // digitised due to noise or whatever.
558  float totalSimADC = 0;
559  for (const auto& iAssociationInfo : associationInfo)
560  totalSimADC += iAssociationInfo.contributionToADC;
561  // Now I know that I can loop again and create the links
562  for (const auto& iAssociationInfo : associationInfo) {
563  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
564  // expected counting to start at 1, not 0. Now changed.
565  outLink.push_back(StripDigiSimLink(channel,
566  iAssociationInfo.trackID,
567  iAssociationInfo.simHitGlobalIndex,
568  iAssociationInfo.tofBin,
569  iAssociationInfo.eventID,
570  iAssociationInfo.contributionToADC / totalSimADC));
571  } // end of loop over associationInfo
572  } // end of loop over the digis
573  } // end of check that iAssociationInfoByChannel is a valid iterator
574 
575  outrawdigi.data = rawdigis;
576 
577  //}
578  }
579 
580  // Now that I've finished with this entry in the map of associations, I can remove it.
581  // Note that there might not be an association if the ADC reading is from noise in which
582  // case associationIsValid will be false.
583  if (iAssociationInfoByChannel != associationInfoForDetId_.end())
584  associationInfoForDetId_.erase(iAssociationInfoByChannel);
585 }
#define LogDebug(id)
unsigned short range
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
const std::unique_ptr< SiPileUpSignals > theSiPileUpSignals
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 > &, std::vector< std::pair< int, std::bitset< 6 >>> &theAffectedAPVvector, CLHEP::HepRandomEngine *)
edm::ESHandle< SiStripLorentzAngle > lorentzAngleHandle
void push_back(const T &t)
Definition: DetSet.h:68
std::map< int, float > mapOfAPVprobabilities
AssociationInfoForDetId associationInfoForDetId_
Structure that holds the information on the SimTrack contributions. Only filled if makeDigiSimLinks_ ...
const std::vector< float > & getMix_TrueInteractions() const
SiDigitalConverter::DigitalRawVecType DigitalRawVecType
SiDigitalConverter::DigitalVecType DigitalVecType
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const std::vector< int > & getMix_bunchCrossing() const
std::vector< unsigned int >::const_iterator ContainerIterator
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
void initializeEvent(const edm::EventSetup &iSetup)
std::pair< ContainerIterator, ContainerIterator > Range
void initializeDetUnit(StripGeomDetUnit const *det, const edm::EventSetup &iSetup)
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:42
const std::unique_ptr< const SiGaussianTailNoiseAdder > theSiNoiseAdder
const std::unique_ptr< SiHitDigitizer > theSiHitDigitizer
std::map< unsigned int, size_t > firstChannelsWithSignal
float getPed(const uint16_t &strip, const Range &range) const
std::map< int, Amplitude > SignalMapType
bool getData(T &iHolder) const
Definition: EventSetup.h:111
static float getNoise(uint16_t strip, const Range &range)
Definition: SiStripNoises.h:74
float getLorentzAngle(const uint32_t &) const
const std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
const std::unique_ptr< SiStripFedZeroSuppression > theSiZeroSuppress
static float getStripGain(const uint16_t &strip, const SiStripApvGain::Range &range)
Definition: SiStripGain.h:73
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::pair< ContainerIterator, ContainerIterator > Range
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
void calculateInstlumiScale(PileupMixingContent *puInfo)
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
RealPedestals
NOTE : turning Noise ON/OFF will make a big change Parameters valid only if Noise = True and ZeroSupp...
SiStripDigitizerAlgorithm(const edm::ParameterSet &conf)
collection_type data
Definition: DetSet.h:80
const Range getRange(const uint32_t detID) const
HLT enums.
const Range getRange(const uint32_t detID) const
std::pair< ContainerIterator, ContainerIterator > Range
std::map< unsigned int, size_t > lastChannelsWithSignal
T get() const
Definition: EventSetup.h:71
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 *)
const int & getMix_bunchSpacing() const
std::string fullPath() const
Definition: FileInPath.cc:163
const Range getRange(const uint32_t &detID) const
std::map< int, std::bitset< 6 > > SiStripTrackerAffectedAPVMap
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:50
bool isValid() const
Definition: ESHandle.h:44
double split
Definition: MVATrainer.cc:139
const ParticleDataTable * pdt
data decode(const unsigned int &value) const
const SiStripApvGain::Range getRange(uint32_t detID) const
Definition: SiStripGain.h:71