CMS 3D CMS Logo

PreMixingSiStripWorker.cc
Go to the documentation of this file.
12 
15 
16 //Data Formats
21 
30 
36 
37 #include "CLHEP/Random/RandFlat.h"
38 
41 
42 #include <map>
43 #include <memory>
44 
46 public:
48  ~PreMixingSiStripWorker() override = default;
49 
50  void initializeEvent(edm::Event const& e, edm::EventSetup const& c) override;
51  void addSignals(edm::Event const& e, edm::EventSetup const& es) override;
52  void addPileups(PileUpEventPrincipal const& pep, edm::EventSetup const& es) override;
53  void put(edm::Event& e, edm::EventSetup const& iSetup, std::vector<PileupSummaryInfo> const& ps, int bs) override;
54 
55 private:
56  void DMinitializeDetUnit(StripGeomDetUnit const* det, const edm::EventSetup& iSetup);
57 
58  // data specifiers
59 
60  edm::InputTag SistripLabelSig_; // name given to collection of SiStrip digis
61  edm::InputTag SiStripPileInputTag_; // InputTag for pileup strips
62  std::string SiStripDigiCollectionDM_; // secondary name to be given to new SiStrip digis
63 
64  edm::InputTag SistripAPVLabelSig_; // where to find vector of dead APVs
67 
73  //edm::ESGetToken<SiStripPedestals, SiStripPedestalsRcd> const pedestalToken_;
76 
77  //
78 
79  typedef float Amplitude;
80  typedef std::pair<uint16_t, Amplitude> RawDigi; // Replacement for SiStripDigi with pulse height instead of integer ADC
81  typedef std::vector<SiStripDigi> OneDetectorMap; // maps by strip ID for later combination - can have duplicate strips
82  typedef std::vector<RawDigi> OneDetectorRawMap; // maps by strip ID for later combination - can have duplicate strips
83  typedef std::map<uint32_t, OneDetectorMap> SiGlobalIndex; // map to all data for each detector ID
84  typedef std::vector<std::pair<uint32_t, OneDetectorRawMap>> SiGlobalRawIndex; // map to all data for each detector ID
85 
87 
90 
91  // variables for temporary storage of mixed hits:
92  typedef std::vector<std::pair<int, Amplitude>> SignalMapType;
93  typedef std::map<uint32_t, SignalMapType> signalMaps;
94 
95  SignalMapType const* getSignal(uint32_t detID) const {
96  auto where = signals_.find(detID);
97  if (where == signals_.end()) {
98  return nullptr;
99  }
100  return &where->second;
101  }
102 
104 
105  // to keep track of dead APVs from HIP interactions
106  typedef std::multimap<uint32_t, std::bitset<6>> APVMap;
107 
109 
110  // for noise adding:
111 
113  bool peakMode;
114  double theThreshold;
118 
119  std::unique_ptr<SiGaussianTailNoiseAdder> theSiNoiseAdder;
120  std::unique_ptr<SiStripFedZeroSuppression> theSiZeroSuppress;
121  std::unique_ptr<SiTrivialDigitalConverter> theSiDigitalConverter;
122 
124 
125  // bad channels for each detector ID
126  std::map<unsigned int, std::vector<bool>> allBadChannels;
127  // channels killed by HIP interactions for each detector ID
128  std::map<unsigned int, std::vector<bool>> allHIPChannels;
129  // first and last channel wit signal for each detector ID
130  std::map<unsigned int, size_t> firstChannelsWithSignal;
131  std::map<unsigned int, size_t> lastChannelsWithSignal;
132 
134  const double fracOfEventsToSimAPV_;
135  const double apv_maxResponse_;
136  const double apv_rate_;
137  const double apv_mVPerQ_;
138  const double apv_fCPerElectron_;
139 
140  //----------------------------
141 
143  public:
144  bool operator()(SiStripDigi i, SiStripDigi j) const { return i.strip() < j.strip(); }
145  };
146 
148  public:
149  bool operator()(RawDigi i, RawDigi j) const { return i.first < j.first; }
150  };
151 };
152 
154  edm::ProducesCollector producesCollector,
156  : pDDToken_(iC.esConsumes(edm::ESInputTag("", ps.getParameter<std::string>("GeometryType")))),
157  deadChannelToken_(iC.esConsumes()),
158  gainToken_(iC.esConsumes(edm::ESInputTag("", ps.getParameter<std::string>("Gain")))),
159  noiseToken_(iC.esConsumes()),
160  thresholdToken_(iC.esConsumes()),
161  SingleStripNoise(ps.getParameter<bool>("SingleStripNoise")),
162  peakMode(ps.getParameter<bool>("APVpeakmode")),
163  theThreshold(ps.getParameter<double>("NoiseSigmaThreshold")),
164  theElectronPerADC(ps.getParameter<double>(peakMode ? "electronPerAdcPeak" : "electronPerAdcDec")),
165  APVSaturationFromHIP_(ps.getParameter<bool>("APVSaturationFromHIP")),
166  theFedAlgo(ps.getParameter<int>("FedAlgorithm_PM")),
167  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo, &iC)),
168  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC, false)), // no premixing
169  includeAPVSimulation_(ps.getParameter<bool>("includeAPVSimulation")),
170  fracOfEventsToSimAPV_(ps.getParameter<double>("fracOfEventsToSimAPV")),
171  apv_maxResponse_(ps.getParameter<double>("apv_maxResponse")),
172  apv_rate_(ps.getParameter<double>("apv_rate")),
173  apv_mVPerQ_(ps.getParameter<double>("apv_mVPerQ")),
174  apv_fCPerElectron_(ps.getParameter<double>("apvfCPerElectron"))
175 
176 {
177  // declare the products to produce
178 
179  SistripLabelSig_ = ps.getParameter<edm::InputTag>("SistripLabelSig");
180  SiStripPileInputTag_ = ps.getParameter<edm::InputTag>("SiStripPileInputTag");
181 
182  SiStripDigiCollectionDM_ = ps.getParameter<std::string>("SiStripDigiCollectionDM");
183  SistripAPVListDM_ = ps.getParameter<std::string>("SiStripAPVListDM");
184 
186  producesCollector.produces<bool>(SiStripDigiCollectionDM_ + "SimulatedAPVDynamicGain");
187 
188  if (APVSaturationFromHIP_) {
189  SistripAPVLabelSig_ = ps.getParameter<edm::InputTag>("SistripAPVLabelSig");
190  SiStripAPVPileInputTag_ = ps.getParameter<edm::InputTag>("SistripAPVPileInputTag");
191  iC.consumes<std::vector<std::pair<int, std::bitset<6>>>>(SistripAPVLabelSig_);
192  }
194 
195  if (includeAPVSimulation_) {
196  tTopoToken_ = iC.esConsumes();
197  apvSimulationParametersToken_ = iC.esConsumes();
198  }
199  // clear local storage for this event
200  SiHitStorage_.clear();
201 
203  if (!rng.isAvailable()) {
204  throw cms::Exception("Psiguration") << "SiStripDigitizer requires the RandomNumberGeneratorService\n"
205  "which is not present in the psiguration file. You must add the service\n"
206  "in the configuration file or remove the modules that require it.";
207  }
208 
209  theSiNoiseAdder = std::make_unique<SiGaussianTailNoiseAdder>(theThreshold);
210 }
211 
213  // initialize individual detectors so we can copy real digitization code:
214 
215  pDD = &iSetup.getData(pDDToken_);
216 
217  for (auto iu = pDD->detUnits().begin(); iu != pDD->detUnits().end(); ++iu) {
218  unsigned int detId = (*iu)->geographicalId().rawId();
219  DetId idet = DetId(detId);
220  unsigned int isub = idet.subdetId();
221  if ((isub == StripSubdetector::TIB) || (isub == StripSubdetector::TID) || (isub == StripSubdetector::TOB) ||
222  (isub == StripSubdetector::TEC)) {
223  auto stripdet = dynamic_cast<StripGeomDetUnit const*>((*iu));
224  assert(stripdet != nullptr);
225  DMinitializeDetUnit(stripdet, iSetup);
226  }
227  }
228 }
229 
231  SiStripBadStrip const& deadChannel = iSetup.getData(deadChannelToken_);
232 
233  unsigned int detId = det->geographicalId().rawId();
234  int numStrips = (det->specificTopology()).nstrips();
235 
236  SiStripBadStrip::Range detBadStripRange = deadChannel.getRange(detId);
237  //storing the bad strip of the the module. the module is not removed but just signal put to 0
238  std::vector<bool>& badChannels = allBadChannels[detId];
239  std::vector<bool>& hipChannels = allHIPChannels[detId];
240  badChannels.clear();
241  badChannels.insert(badChannels.begin(), numStrips, false);
242  hipChannels.clear();
243  hipChannels.insert(hipChannels.begin(), numStrips, false);
244 
245  for (SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
246  SiStripBadStrip::data fs = deadChannel.decode(*it);
247  for (int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip)
248  badChannels[strip] = true;
249  }
250  firstChannelsWithSignal[detId] = numStrips;
252 }
253 
255  // fill in maps of hits
256 
258 
259  if (e.getByLabel(SistripLabelSig_, input)) {
260  OneDetectorMap LocalMap;
261 
262  //loop on all detsets (detectorIDs) inside the input collection
264  for (; DSViter != input->end(); DSViter++) {
265 #ifdef DEBUG
266  LogDebug("PreMixingSiStripWorker") << "Processing DetID " << DSViter->id;
267 #endif
268 
269  LocalMap.clear();
270  LocalMap.reserve((DSViter->data).size());
271  LocalMap.insert(LocalMap.end(), (DSViter->data).begin(), (DSViter->data).end());
272 
273  SiHitStorage_.insert(SiGlobalIndex::value_type(DSViter->id, LocalMap));
274  }
275  }
276 
277  // keep here for future reference. In current implementation, HIP killing is done once in PU file
278  /* if(APVSaturationFromHIP_) {
279  edm::Handle<std::vector<std::pair<int,std::bitset<6>> > > APVinput;
280 
281  if( e.getByLabel(SistripAPVLabelSig_,APVinput) ) {
282 
283  std::vector<std::pair<int,std::bitset<6>> >::const_iterator entry = APVinput->begin();
284  for( ; entry != APVinput->end(); entry++) {
285  theAffectedAPVmap_.insert(APVMap::value_type(entry->first, entry->second));
286  }
287  }
288  } */
289 
290 } // end of addSiStripSignals
291 
293  LogDebug("PreMixingSiStripWorker") << "\n===============> adding pileups from event " << pep.principal().id()
294  << " for bunchcrossing " << pep.bunchCrossing();
295 
296  // fill in maps of hits; same code as addSignals, except now applied to the pileup events
297 
299  pep.getByLabel(SiStripPileInputTag_, inputHandle);
300 
301  if (inputHandle.isValid()) {
302  const auto& input = *inputHandle;
303 
304  OneDetectorMap LocalMap;
305 
306  //loop on all detsets (detectorIDs) inside the input collection
308  for (; DSViter != input.end(); DSViter++) {
309 #ifdef DEBUG
310  LogDebug("PreMixingSiStripWorker") << "Pileups: Processing DetID " << DSViter->id;
311 #endif
312 
313  // find correct local map (or new one) for this detector ID
314 
315  SiGlobalIndex::const_iterator itest;
316 
317  itest = SiHitStorage_.find(DSViter->id);
318 
319  if (itest != SiHitStorage_.end()) { // this detID already has hits, add to existing map
320 
321  LocalMap = itest->second;
322 
323  // fill in local map with extra channels
324  LocalMap.insert(LocalMap.end(), (DSViter->data).begin(), (DSViter->data).end());
325  std::stable_sort(LocalMap.begin(), LocalMap.end(), PreMixingSiStripWorker::StrictWeakOrdering());
326  SiHitStorage_[DSViter->id] = LocalMap;
327 
328  } else { // fill local storage with this information, put in global collection
329 
330  LocalMap.clear();
331  LocalMap.reserve((DSViter->data).size());
332  LocalMap.insert(LocalMap.end(), (DSViter->data).begin(), (DSViter->data).end());
333 
334  SiHitStorage_.insert(SiGlobalIndex::value_type(DSViter->id, LocalMap));
335  }
336  }
337 
338  if (APVSaturationFromHIP_) {
340  pep.getByLabel(SiStripAPVPileInputTag_, inputAPVHandle);
341 
342  if (inputAPVHandle.isValid()) {
343  const auto& APVinput = inputAPVHandle;
344 
345  std::vector<std::pair<int, std::bitset<6>>>::const_iterator entry = APVinput->begin();
346  for (; entry != APVinput->end(); entry++) {
347  theAffectedAPVmap_.insert(APVMap::value_type(entry->first, entry->second));
348  }
349  }
350  }
351  }
352 }
353 
355  edm::EventSetup const& iSetup,
356  std::vector<PileupSummaryInfo> const& ps,
357  int bs) {
358  // set up machinery to do proper noise adding:
359  SiStripGain const& gain = iSetup.getData(gainToken_);
360  SiStripNoises const& noise = iSetup.getData(noiseToken_);
362  //SiStripPedestals const& pedestal = iSetup.getData(pedestalToken_);
363  SiStripApvSimulationParameters const* apvSimulationParameters = nullptr;
364  TrackerTopology const* tTopo = nullptr;
365 
367  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
368 
369  const bool simulateAPVInThisEvent = includeAPVSimulation_ && (CLHEP::RandFlat::shoot(engine) < fracOfEventsToSimAPV_);
370  float nTruePU = 0.; // = ps.getTrueNumInteractions();
371  if (simulateAPVInThisEvent) {
372  tTopo = &iSetup.getData(tTopoToken_);
373  apvSimulationParameters = &iSetup.getData(apvSimulationParametersToken_);
374  const auto it = std::find_if(
375  std::begin(ps), std::end(ps), [](const PileupSummaryInfo& bxps) { return bxps.getBunchCrossing() == 0; });
376  if (it != std::begin(ps)) {
377  nTruePU = it->getTrueNumInteractions();
378  } else {
379  edm::LogWarning("PreMixingSiStripWorker") << "Could not find PileupSummaryInfo for current bunch crossing";
380  nTruePU = std::begin(ps)->getTrueNumInteractions();
381  }
382  }
383 
384  std::map<int, std::bitset<6>> DeadAPVList;
385 
386  // First, have to convert all ADC counts to raw pulse heights so that values can be added properly
387  // In PreMixing, pulse heights are saved with ADC = sqrt(9.0*PulseHeight) - have to undo.
388 
389  // This is done here because it's the only place we have access to EventSetup
390  // Simultaneously, merge lists of hit channels in each DetId.
391  // Signal Digis are in the list first, have to merge lists of hit strips on the fly,
392  // add signals on duplicates later
393 
394  // Now, loop over hits and add them to the map in the proper sorted order
395  // Note: We are assuming that the hits from the Signal events have been created in
396  // "PreMix" mode, rather than in the standard ADC conversion routines. If not, this
397  // doesn't work at all.
398 
399  // At the moment, both Signal and Reconstituted PU hits have the same compression algorithm.
400  // If this were different, and one needed gains, the conversion back to pulse height can only
401  // be done in this routine. So, yes, there is an extra loop over hits here in the current code,
402  // because, in principle, one could convert to pulse height during the read/store phase.
403 
404  for (SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin(); IDet != SiHitStorage_.end(); IDet++) {
405  uint32_t detID = IDet->first;
406 
407  auto const& LocalMap = IDet->second;
408 
409  //loop over hit strips for this DetId, do conversion to pulse height, store.
410  SiRawDigis_.emplace_back(detID, OneDetectorRawMap());
411  auto& localRawMap = SiRawDigis_.back().second;
412  localRawMap.reserve(LocalMap.size());
413  for (auto const& iLocal : LocalMap) {
414  uint16_t currentStrip = iLocal.strip();
415  float signal = float(iLocal.adc());
416  if (iLocal.adc() == 1022)
417  signal = 1500.f; // average values for overflows
418  if (iLocal.adc() == 1023)
419  signal = 3000.f;
420 
421  //convert signals back to raw counts
422  constexpr float inv9 = 1.f / 9.0f;
423  float ReSignal = signal * signal * inv9; // The PreMixing conversion is adc = sqrt(9.0*pulseHeight)
424 
425  // RawDigi NewRawDigi = std::make_pair(currentStrip, ReSignal);
426 
427  localRawMap.emplace_back(currentStrip, ReSignal);
428  }
429  }
430 
431  // If we are killing APVs, merge list of dead ones before we digitize
432 
433  int NumberOfBxBetweenHIPandEvent = 1e3;
434 
435  if (APVSaturationFromHIP_) {
436  // calculate affected BX parameter
437 
438  bool HasAtleastOneAffectedAPV = false;
439  while (!HasAtleastOneAffectedAPV) {
440  for (int bx = floor(300.0 / 25.0); bx > 0; bx--) { //Reminder: make these numbers not hard coded!!
441  float temp = CLHEP::RandFlat::shoot(engine) < 0.5 ? 1 : 0;
442  if (temp == 1 && bx < NumberOfBxBetweenHIPandEvent) {
443  NumberOfBxBetweenHIPandEvent = bx;
444  HasAtleastOneAffectedAPV = true;
445  }
446  }
447  }
448 
449  APVMap::const_iterator iAPVchk;
450  uint32_t formerID = 0;
451  uint32_t currentID;
452  std::bitset<6> NewAPVBits;
453 
454  for (APVMap::const_iterator iAPV = theAffectedAPVmap_.begin(); iAPV != theAffectedAPVmap_.end(); ++iAPV) {
455  currentID = iAPV->first;
456 
457  if (currentID == formerID) { // we have to OR these
458  for (int ibit = 0; ibit < 6; ++ibit) {
459  NewAPVBits[ibit] = NewAPVBits[ibit] || (iAPV->second)[ibit];
460  }
461  } else {
462  DeadAPVList[currentID] = NewAPVBits;
463  //save pointers for next iteration
464  formerID = currentID;
465  NewAPVBits = iAPV->second;
466  }
467 
468  iAPVchk = iAPV;
469  if ((++iAPVchk) == theAffectedAPVmap_.end()) { //make sure not to lose the last one
470  DeadAPVList[currentID] = NewAPVBits;
471  }
472  }
473  }
474  //
475 
476  // Ok, done with merging raw signals and APVs - now add signals on duplicate strips
477 
478  // collection of Digis to put in the event
479  std::vector<edm::DetSet<SiStripDigi>> vSiStripDigi;
480 
481  // loop through our collection of detectors, merging hits and making a new list of "signal" digis
482 
483  // clear some temporary storage for later digitization:
484 
485  signals_.clear();
486 
487  // big loop over Detector IDs:
488  for (auto const& IDet : SiRawDigis_) {
489  uint32_t detID = IDet.first;
490 
491  // create merged map:
492  auto& lSignals = signals_[detID];
493 
494  auto const& LocalMap = IDet.second;
495 
496  // guess max occupancy...
497  lSignals.reserve(std::min(LocalMap.size(), 512ul));
498 
499  //counter variables
500  int formerStrip = -1;
501  int currentStrip;
502  float ADCSum = 0;
503 
504  //loop over hit strips for this DetId, add duplicates
505 
506  OneDetectorRawMap::const_iterator iLocalchk;
507  OneDetectorRawMap::const_iterator iLocal = LocalMap.begin();
508  for (; iLocal != LocalMap.end(); ++iLocal) {
509  currentStrip = iLocal->first; // strip is first element
510 
511  if (currentStrip == formerStrip) { // we have to add these digis together
512 
513  ADCSum += iLocal->second; // raw pulse height is second element.
514  } else {
515  if (formerStrip != -1) {
516  lSignals.emplace_back(formerStrip, ADCSum);
517  }
518  // save pointers for next iteration
519  formerStrip = currentStrip;
520  ADCSum = iLocal->second; // lone ADC
521  }
522 
523  iLocalchk = iLocal;
524  if ((++iLocalchk) == LocalMap.end()) { //make sure not to lose the last one
525  lSignals.emplace_back(formerStrip, ADCSum);
526  }
527  }
528  }
529 
530  //Now, do noise, zero suppression, take into account bad channels, etc.
531  // This section stolen from SiStripDigitizerAlgorithm
532  // must loop over all detIds in the tracker to get all of the noise added properly.
533  for (const auto& iu : pDD->detUnits()) {
534  const StripGeomDetUnit* sgd = dynamic_cast<const StripGeomDetUnit*>(iu);
535  if (sgd != nullptr) {
536  uint32_t detID = sgd->geographicalId().rawId();
537 
538  edm::DetSet<SiStripDigi> SSD(detID); // Make empty collection with this detector ID
539 
540  int numStrips = (sgd->specificTopology()).nstrips();
541 
542  // see if there is some signal on this detector
543 
544  auto const* lSignal = getSignal(detID);
545 
546  std::vector<float> detAmpl(numStrips, 0.);
547  if (lSignal) {
548  for (const auto& amp : *lSignal) {
549  detAmpl[amp.first] = amp.second;
550  }
551  }
552 
553  //removing signal from the dead (and HIP effected) strips
554  std::vector<bool>& badChannels = allBadChannels[detID];
555 
556  for (int strip = 0; strip < numStrips; ++strip) {
557  if (badChannels[strip])
558  detAmpl[strip] = 0.;
559  }
560 
561  if (simulateAPVInThisEvent) {
562  // Get index in apv baseline distributions corresponding to z of detSet and PU
563  const StripTopology* topol = dynamic_cast<const StripTopology*>(&(sgd->specificTopology()));
564  LocalPoint localPos = topol->localPosition(0);
565  GlobalPoint globalPos = sgd->surface().toGlobal(Local3DPoint(localPos.x(), localPos.y(), localPos.z()));
566  float detSet_z = fabs(globalPos.z());
567  float detSet_r = globalPos.perp();
568 
569  const uint32_t SubDet = DetId(detID).subdetId();
570  // Simulate APV response for each strip
571  for (int strip = 0; strip < numStrips; ++strip) {
572  if (detAmpl[strip] > 0) {
573  // Convert charge from electrons to fC
574  double stripCharge = detAmpl[strip] * apv_fCPerElectron_;
575 
576  // Get APV baseline
577  double baselineV = 0;
579  baselineV = apvSimulationParameters->sampleTIB(tTopo->tibLayer(detID), detSet_z, nTruePU, engine);
580  } else if (SubDet == SiStripSubdetector::TOB) {
581  baselineV = apvSimulationParameters->sampleTOB(tTopo->tobLayer(detID), detSet_z, nTruePU, engine);
582  } else if (SubDet == SiStripSubdetector::TID) {
583  baselineV = apvSimulationParameters->sampleTID(tTopo->tidWheel(detID), detSet_r, nTruePU, engine);
584  } else if (SubDet == SiStripSubdetector::TEC) {
585  baselineV = apvSimulationParameters->sampleTEC(tTopo->tecWheel(detID), detSet_r, nTruePU, engine);
586  }
587  // Fitted parameters from G Hall/M Raymond
588  double maxResponse = apv_maxResponse_;
589  double rate = apv_rate_;
590 
591  double outputChargeInADC = 0;
592  if (baselineV < apv_maxResponse_) {
593  // Convert V0 into baseline charge
594  double baselineQ = -1.0 * rate * log(2 * maxResponse / (baselineV + maxResponse) - 1);
595 
596  // Add charge deposited in this BX
597  double newStripCharge = baselineQ + stripCharge;
598 
599  // Apply APV response
600  double signalV = 2 * maxResponse / (1 + exp(-1.0 * newStripCharge / rate)) - maxResponse;
601  double gain = signalV - baselineV;
602 
603  // Convert gain (mV) to charge (assuming linear region of APV) and then to electrons
604  double outputCharge = gain / apv_mVPerQ_;
605  outputChargeInADC = outputCharge / apv_fCPerElectron_;
606  }
607 
608  // Output charge back to original container
609  detAmpl[strip] = outputChargeInADC;
610  }
611  }
612  }
613 
614  if (APVSaturationFromHIP_) {
615  std::bitset<6>& bs = DeadAPVList[detID];
616 
617  if (bs.any()) {
618  // Here below is the scaling function which describes the evolution of the baseline (i.e. how the charge is suppressed).
619  // This must be replaced as soon as we have a proper modeling of the baseline evolution from VR runs
620  float Shift =
621  1 - NumberOfBxBetweenHIPandEvent / floor(300.0 / 25.0); //Reminder: make these numbers not hardcoded!!
622  float randomX = CLHEP::RandFlat::shoot(engine);
623  float scalingValue = (randomX - Shift) * 10.0 / 7.0 - 3.0 / 7.0;
624 
625  for (int strip = 0; strip < numStrips; ++strip) {
626  if (!badChannels[strip] && bs[strip / 128] == 1) {
627  detAmpl[strip] *= scalingValue > 0 ? scalingValue : 0.0;
628  }
629  }
630  }
631  }
632 
633  SiStripNoises::Range detNoiseRange = noise.getRange(detID);
634  SiStripApvGain::Range detGainRange = gain.getRange(detID);
635 
636  // Gain conversion is already done during signal adding
637  //convert our signals back to raw counts so that we can add noise properly:
638 
639  /*
640  if(theSignal) {
641  for(unsigned int iv = 0; iv!=detAmpl.size(); iv++) {
642  float signal = detAmpl[iv];
643  if(signal > 0) {
644  float gainValue = gain.getStripGain(iv, detGainRange);
645  signal *= theElectronPerADC/gainValue;
646  detAmpl[iv] = signal;
647  }
648  }
649  }
650  */
651 
652  //SiStripPedestals::Range detPedestalRange = pedestal.getRange(detID);
653 
654  // -----------------------------------------------------------
655 
656  size_t firstChannelWithSignal = 0;
657  size_t lastChannelWithSignal = numStrips;
658 
659  if (SingleStripNoise) {
660  std::vector<float> noiseRMSv;
661  noiseRMSv.clear();
662  noiseRMSv.insert(noiseRMSv.begin(), numStrips, 0.);
663  for (int strip = 0; strip < numStrips; ++strip) {
664  if (!badChannels[strip]) {
665  float gainValue = gain.getStripGain(strip, detGainRange);
666  noiseRMSv[strip] = (noise.getNoise(strip, detNoiseRange)) * theElectronPerADC / gainValue;
667  }
668  }
669  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
670  } else {
671  int RefStrip = int(numStrips / 2.);
672  while (RefStrip < numStrips &&
673  badChannels[RefStrip]) { //if the refstrip is bad, I move up to when I don't find it
674  RefStrip++;
675  }
676  if (RefStrip < numStrips) {
677  float RefgainValue = gain.getStripGain(RefStrip, detGainRange);
678  float RefnoiseRMS = noise.getNoise(RefStrip, detNoiseRange) * theElectronPerADC / RefgainValue;
679 
680  theSiNoiseAdder->addNoise(
681  detAmpl, firstChannelWithSignal, lastChannelWithSignal, numStrips, RefnoiseRMS, engine);
682  }
683  }
684 
685  auto const& data = theSiDigitalConverter->convert(detAmpl, &gain, detID);
686  // if suppression is trival just copy data
687  if (5 == theFedAlgo)
688  SSD.data = data;
689  else
690  theSiZeroSuppress->suppress(data, SSD.data, detID, noise, threshold);
691 
692  LogDebug("PreMixingSiStripWorker") << detID << " " << data.size() << ' ' << SSD.data.size();
693 
694  // stick this into the global vector of detector info
695  vSiStripDigi.push_back(SSD);
696 
697  } // end of loop over one detector
698 
699  } // end of big loop over all detector IDs
700 
701  // put the collection of digis in the event
702  edm::LogInfo("PreMixingSiStripWorker") << "total # Merged strips: " << vSiStripDigi.size();
703 
704  // make new digi collection
705 
706  std::unique_ptr<edm::DetSetVector<SiStripDigi>> MySiStripDigis(new edm::DetSetVector<SiStripDigi>(vSiStripDigi));
707 
708  // put collection
709 
710  e.put(std::move(MySiStripDigis), SiStripDigiCollectionDM_);
711  e.put(std::make_unique<bool>(simulateAPVInThisEvent), SiStripDigiCollectionDM_ + "SimulatedAPVDynamicGain");
712 
713  // clear local storage for this event
714  SiHitStorage_.clear();
715  SiRawDigis_.clear();
716  signals_.clear();
717 }
718 
bool operator()(SiStripDigi i, SiStripDigi j) const
std::map< unsigned int, std::vector< bool > > allBadChannels
static constexpr auto TEC
void addSignals(edm::Event const &e, edm::EventSetup const &es) override
std::vector< SiStripDigi > DigitalVecType
std::map< unsigned int, std::vector< bool > > allHIPChannels
PreMixingSiStripWorker(const edm::ParameterSet &ps, edm::ProducesCollector, edm::ConsumesCollector &&iC)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned int tobLayer(const DetId &id) const
T perp() const
Definition: PV3DBase.h:69
float sampleTOB(layerid layer, float z, float pu, CLHEP::HepRandomEngine *engine) const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
ProductRegistryHelper::BranchAliasSetterT< ProductType > produces()
std::map< unsigned int, size_t > lastChannelsWithSignal
T z() const
Definition: PV3DBase.h:61
void reserve(size_t s)
Definition: DetSetVector.h:145
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > const pDDToken_
std::vector< unsigned int >::const_iterator ContainerIterator
unsigned int tidWheel(const DetId &id) const
unsigned int tecWheel(const DetId &id) const
std::map< uint32_t, OneDetectorMap > SiGlobalIndex
void DMinitializeDetUnit(StripGeomDetUnit const *det, const edm::EventSetup &iSetup)
float sampleTIB(layerid layer, float z, float pu, CLHEP::HepRandomEngine *engine) const
edm::ESGetToken< SiStripNoises, SiStripNoisesRcd > const noiseToken_
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
assert(be >=bs)
const Range getRange(const uint32_t detID) const
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
std::multimap< uint32_t, std::bitset< 6 > > APVMap
static std::string const input
Definition: EdmProvDump.cc:50
edm::EventPrincipal const & principal()
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
std::vector< SiStripDigi > OneDetectorMap
std::unique_ptr< SiGaussianTailNoiseAdder > theSiNoiseAdder
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
std::pair< ContainerIterator, ContainerIterator > Range
double f[11][100]
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
~PreMixingSiStripWorker() override=default
std::vector< std::pair< uint32_t, OneDetectorRawMap > > SiGlobalRawIndex
A Digi for the silicon strip detector, containing both strip and adc information, and suitable for st...
Definition: SiStripDigi.h:12
static constexpr auto TOB
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
std::unique_ptr< SiStripFedZeroSuppression > theSiZeroSuppress
void put(edm::Event &e, edm::EventSetup const &iSetup, std::vector< PileupSummaryInfo > const &ps, int bs) override
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
Log< level::Info, false > LogInfo
edm::ESGetToken< SiStripBadStrip, SiStripBadChannelRcd > const deadChannelToken_
Definition: DetId.h:17
const TrackerGeometry * pDD
float sampleTID(layerid wheel, float r, float pu, CLHEP::HepRandomEngine *engine) const
static constexpr auto TIB
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
void initializeEvent(edm::Event const &e, edm::EventSetup const &c) override
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const int getBunchCrossing() const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
SignalMapType const * getSignal(uint32_t detID) const
std::pair< uint16_t, Amplitude > RawDigi
bool isValid() const
Definition: HandleBase.h:70
double rate(double x)
Definition: Constants.cc:3
collection_type data
Definition: DetSet.h:80
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
std::pair< ContainerIterator, ContainerIterator > Range
std::map< unsigned int, size_t > firstChannelsWithSignal
data decode(const unsigned int &value) const
virtual LocalPoint localPosition(float strip) const =0
void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &es) override
bool isAvailable() const
Definition: Service.h:40
unsigned int tibLayer(const DetId &id) const
std::vector< RawDigi > OneDetectorRawMap
Log< level::Warning, false > LogWarning
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:47
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
EventID const & id() const
float sampleTEC(layerid wheel, float r, float pu, CLHEP::HepRandomEngine *engine) const
static constexpr auto TID
#define DEFINE_PREMIXING_WORKER(TYPE)
std::map< uint32_t, SignalMapType > signalMaps
edm::ESGetToken< SiStripGain, SiStripGainSimRcd > const gainToken_
def move(src, dest)
Definition: eostools.py:511
edm::ESGetToken< SiStripThreshold, SiStripThresholdRcd > const thresholdToken_
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
std::vector< std::pair< int, Amplitude > > SignalMapType
edm::ESGetToken< SiStripApvSimulationParameters, SiStripApvSimulationParametersRcd > apvSimulationParametersToken_
#define LogDebug(id)
SiDigitalConverter::DigitalVecType DigitalVecType
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9