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