CMS 3D CMS Logo

SiStripAPVRestorer.cc
Go to the documentation of this file.
4 
5 #include <cmath>
6 #include <iostream>
7 #include <algorithm>
8 
10  : quality_cache_id(-1),
11  noise_cache_id(-1),
12  pedestal_cache_id(-1),
13  forceNoRestore_(conf.getParameter<bool>("ForceNoRestore")),
14  inspectAlgo_(conf.getParameter<std::string>("APVInspectMode")),
15  restoreAlgo_(conf.getParameter<std::string>("APVRestoreMode")),
16  // CM map settings
17  useRealMeanCM_(conf.getParameter<bool>("useRealMeanCM")),
18  meanCM_(conf.getParameter<int32_t>("MeanCM")),
19  deltaCMThreshold_(conf.getParameter<uint32_t>("DeltaCMThreshold")),
20  // abnormal baseline inspect
21  fraction_(conf.getParameter<double>("Fraction")),
22  deviation_(conf.getParameter<uint32_t>("Deviation")),
23  // nullInspect
24  restoreThreshold_(conf.getParameter<double>("restoreThreshold")),
25  // baseline and saturation inspect
26  nSaturatedStrip_(conf.getParameter<uint32_t>("nSaturatedStrip")),
27  // FlatRegionsFinder
28  nSigmaNoiseDerTh_(conf.getParameter<uint32_t>("nSigmaNoiseDerTh")),
29  consecThreshold_(conf.getParameter<uint32_t>("consecThreshold")),
30  nSmooth_(conf.getParameter<uint32_t>("nSmooth")),
31  distortionThreshold_(conf.getParameter<uint32_t>("distortionThreshold")),
32  applyBaselineCleaner_(conf.getParameter<bool>("ApplyBaselineCleaner")),
33  cleaningSequence_(conf.getParameter<uint32_t>("CleaningSequence")),
34  // cleaner_localmin
35  slopeX_(conf.getParameter<int32_t>("slopeX")),
36  slopeY_(conf.getParameter<int32_t>("slopeY")),
37  // cleaner_monotony
38  hitStripThreshold_(conf.getParameter<uint32_t>("hitStripThreshold")),
39  // baseline follower (restore)
40  minStripsToFit_(conf.getParameter<uint32_t>("minStripsToFit")),
41  applyBaselineRejection_(conf.getParameter<bool>("ApplyBaselineRejection")),
42  filteredBaselineMax_(conf.getParameter<double>("filteredBaselineMax")),
43  filteredBaselineDerivativeSumSquare_(conf.getParameter<double>("filteredBaselineDerivativeSumSquare")),
44  // derivative follower restorer
45  gradient_threshold_(conf.getParameter<int>("discontinuityThreshold")),
46  last_gradient_(conf.getParameter<int>("lastGradient")),
47  size_window_(conf.getParameter<int>("sizeWindow")),
48  width_cluster_(conf.getParameter<int>("widthCluster")) {
49  if (restoreAlgo_ == "BaselineFollower" && inspectAlgo_ != "BaselineFollower" && inspectAlgo_ != "Hybrid")
50  throw cms::Exception("Incompatible Algorithm")
51  << "The BaselineFollower restore method requires the BaselineFollower (or Hybrid) inspect method";
52 }
53 
55  uint32_t n_cache_id = es.get<SiStripNoisesRcd>().cacheIdentifier();
56  uint32_t q_cache_id = es.get<SiStripQualityRcd>().cacheIdentifier();
57  uint32_t p_cache_id = es.get<SiStripPedestalsRcd>().cacheIdentifier();
58 
59  if (n_cache_id != noise_cache_id) {
61  noise_cache_id = n_cache_id;
62  } else {
63  noise_cache_id = n_cache_id;
64  }
65  if (q_cache_id != quality_cache_id) {
67  quality_cache_id = q_cache_id;
68  } else {
69  quality_cache_id = q_cache_id;
70  }
71 
72  if (p_cache_id != pedestal_cache_id) {
74  pedestal_cache_id = p_cache_id;
75  } else {
76  pedestal_cache_id = p_cache_id;
77  }
78 }
79 
80 uint16_t SiStripAPVRestorer::inspectAndRestore(uint32_t detId,
81  uint16_t firstAPV,
82  const digivector_t& rawDigisPedSubtracted,
83  digivector_t& processedRawDigi,
84  const medians_t& vmedians) {
85  auto nAPVFlagged = inspect(detId, firstAPV, rawDigisPedSubtracted, vmedians);
86  restore(firstAPV, processedRawDigi);
87  return nAPVFlagged;
88 }
89 
90 uint16_t SiStripAPVRestorer::inspect(uint32_t detId,
91  uint16_t firstAPV,
92  const digivector_t& digis,
93  const medians_t& vmedians) {
94  detId_ = detId;
95 
96  apvFlags_.assign(6, "");
97  apvFlagsBoolOverride_.assign(6, false);
98  apvFlagsBool_.clear();
99  median_.assign(6, -999);
100  badAPVs_.assign(6, false);
101  smoothedMaps_.clear();
102  baselineMap_.clear();
103  for (const auto& med : vmedians) {
104  auto iAPV = med.first;
105  median_[iAPV] = med.second;
106  badAPVs_[iAPV] = qualityHandle->IsApvBad(detId_, iAPV);
107  }
108 
109  uint16_t nAPVFlagged = 0;
110  if (inspectAlgo_ == "Hybrid")
111  nAPVFlagged = hybridFormatInspect(firstAPV, digis);
112  else if (inspectAlgo_ == "HybridEmulation")
113  nAPVFlagged = hybridEmulationInspect(firstAPV, digis);
114  else if (inspectAlgo_ == "BaselineFollower")
115  nAPVFlagged = baselineFollowerInspect(firstAPV, digis);
116  else if (inspectAlgo_ == "AbnormalBaseline")
117  nAPVFlagged = abnormalBaselineInspect(firstAPV, digis);
118  else if (inspectAlgo_ == "Null")
119  nAPVFlagged = nullInspect(firstAPV, digis);
120  else if (inspectAlgo_ == "BaselineAndSaturation")
121  nAPVFlagged = baselineAndSaturationInspect(firstAPV, digis);
122  else if (inspectAlgo_ == "DerivativeFollower")
123  nAPVFlagged = forceRestoreInspect(firstAPV, digis);
124  else
125  throw cms::Exception("Unregistered Inspect Algorithm")
126  << "SiStripAPVRestorer possibilities: (Null), (AbnormalBaseline), (BaselineFollower), (BaselineAndSaturation), "
127  "(DerivativeFollower), (Hybrid), (HybridEmulation)";
128 
129  // update boolean version of APV flags
130  for (size_t i = 0; i < apvFlags_.size(); ++i) {
131  apvFlagsBool_.push_back((!apvFlags_[i].empty()) && (!apvFlagsBoolOverride_[i]));
132  }
133 
134  return nAPVFlagged;
135 }
136 
137 void SiStripAPVRestorer::restore(uint16_t firstAPV, digivector_t& digis) {
138  if (forceNoRestore_)
139  return;
140 
141  for (uint16_t iAPV = firstAPV; iAPV < digis.size() / nTotStripsPerAPV + firstAPV; ++iAPV) {
142  auto algoToUse = apvFlags_[iAPV];
143  if (!algoToUse.empty()) {
144  if (algoToUse == "Flat") {
145  flatRestore(iAPV, firstAPV, digis);
146  } else if (algoToUse == "BaselineFollower") {
147  baselineFollowerRestore(iAPV, firstAPV, median_[iAPV], digis);
148  } else if (algoToUse == "DerivativeFollower") {
149  derivativeFollowerRestore(iAPV, firstAPV, digis);
150  } else {
151  throw cms::Exception("Unregistered Restore Algorithm")
152  << "SiStripAPVRestorer possibilities: (Flat), (BaselineFollower), (DerivativeFollower)";
153  }
154  }
155  }
156 }
157 
158 // Inspect method implementations
159 
160 inline uint16_t SiStripAPVRestorer::hybridFormatInspect(uint16_t firstAPV, const digivector_t& digis) {
161  digivector_t singleAPVdigi;
162  singleAPVdigi.reserve(nTotStripsPerAPV);
163  uint16_t nAPVflagged = 0;
164  for (uint16_t iAPV = firstAPV; iAPV < digis.size() / nTotStripsPerAPV + firstAPV; ++iAPV) {
165  digimap_t smoothedmap;
166  if (!badAPVs_[iAPV]) {
167  uint16_t stripsPerAPV = 0;
168  singleAPVdigi.clear();
169  for (int16_t strip = (iAPV - firstAPV) * nTotStripsPerAPV; strip < (iAPV - firstAPV + 1) * nTotStripsPerAPV;
170  ++strip) {
171  const uint16_t adc = digis[strip];
172  singleAPVdigi.push_back(adc);
173  if (adc > 0)
174  ++stripsPerAPV;
175  }
176 
177  if (stripsPerAPV > 64) {
178  flatRegionsFinder(singleAPVdigi, smoothedmap, iAPV);
179  apvFlags_[iAPV] = restoreAlgo_; //specify any algo to make the restore
180  nAPVflagged++;
181  }
182  }
183  smoothedMaps_.emplace_hint(smoothedMaps_.end(), iAPV, std::move(smoothedmap));
184  }
185  return nAPVflagged;
186 }
187 
188 inline uint16_t SiStripAPVRestorer::baselineFollowerInspect(uint16_t firstAPV, const digivector_t& digis) {
189  digivector_t singleAPVdigi;
190  singleAPVdigi.reserve(nTotStripsPerAPV);
191  uint16_t nAPVflagged = 0;
192 
193  auto itCMMap = std::end(meanCMmap_);
194  if (useRealMeanCM_)
195  itCMMap = meanCMmap_.find(detId_);
196 
197  for (uint16_t iAPV = firstAPV; iAPV < firstAPV + digis.size() / nTotStripsPerAPV; ++iAPV) {
198  digimap_t smoothedmap;
199  if (!badAPVs_[iAPV]) {
200  float MeanAPVCM = meanCM_;
201  if (useRealMeanCM_ && (std::end(meanCMmap_) != itCMMap))
202  MeanAPVCM = itCMMap->second[iAPV];
203 
204  singleAPVdigi.clear();
205  std::copy_n(std::begin(digis) + nTotStripsPerAPV * (iAPV - firstAPV),
207  std::back_inserter(singleAPVdigi));
208 
209  const float DeltaCM = median_[iAPV] - MeanAPVCM;
210  if ((DeltaCM < 0) && (std::abs(DeltaCM) > deltaCMThreshold_)) {
211  if (!flatRegionsFinder(singleAPVdigi, smoothedmap, iAPV)) {
212  apvFlags_[iAPV] = restoreAlgo_; //specify any algo to make the restore
213  ++nAPVflagged;
214  }
215  }
216  }
217  smoothedMaps_.emplace_hint(smoothedMaps_.end(), iAPV, std::move(smoothedmap));
218  }
219  return nAPVflagged;
220 }
221 
222 inline uint16_t SiStripAPVRestorer::forceRestoreInspect(uint16_t firstAPV, const digivector_t& digis) {
223  uint16_t nAPVflagged = 0;
224  for (uint16_t iAPV = firstAPV; iAPV < firstAPV + digis.size() / nTotStripsPerAPV; ++iAPV) {
225  if (!badAPVs_[iAPV]) {
226  apvFlags_[iAPV] = restoreAlgo_; //specify any algo to make the restore
227  nAPVflagged++;
228  }
229  }
230  return nAPVflagged;
231 }
232 
233 inline uint16_t SiStripAPVRestorer::baselineAndSaturationInspect(uint16_t firstAPV, const digivector_t& digis) {
234  digivector_t singleAPVdigi;
235  singleAPVdigi.reserve(nTotStripsPerAPV);
236  uint16_t nAPVflagged = 0;
237 
238  auto itCMMap = std::end(meanCMmap_);
239  if (useRealMeanCM_)
240  itCMMap = meanCMmap_.find(detId_);
241 
242  for (uint16_t iAPV = firstAPV; iAPV < firstAPV + digis.size() / nTotStripsPerAPV; ++iAPV) {
243  if (!badAPVs_[iAPV]) {
244  float MeanAPVCM = meanCM_;
245  if (useRealMeanCM_ && (std::end(meanCMmap_) != itCMMap))
246  MeanAPVCM = itCMMap->second[iAPV];
247 
248  singleAPVdigi.clear();
249  uint16_t nSatStrip = 0;
250  for (int16_t strip = (iAPV - firstAPV) * nTotStripsPerAPV; strip < (iAPV - firstAPV + 1) * nTotStripsPerAPV;
251  ++strip) {
252  const uint16_t digi = digis[strip];
253  singleAPVdigi.push_back(digi);
254  if (digi >= 1023)
255  ++nSatStrip;
256  }
257 
258  const float DeltaCM = median_[iAPV] - MeanAPVCM;
259  if ((DeltaCM < 0) && (std::abs(DeltaCM) > deltaCMThreshold_) && (nSatStrip >= nSaturatedStrip_)) {
260  apvFlags_[iAPV] = restoreAlgo_; //specify any algo to make the restore
261  ++nAPVflagged;
262  }
263  }
264  }
265  return nAPVflagged;
266 }
267 
268 inline uint16_t SiStripAPVRestorer::abnormalBaselineInspect(uint16_t firstAPV, const digivector_t& digis) {
270 
271  uint16_t nAPVflagged = 0;
272 
273  auto itCMMap = std::end(meanCMmap_);
274  if (useRealMeanCM_)
275  itCMMap = meanCMmap_.find(detId_);
276 
277  int devCount = 0, qualityCount = 0, minstrip = 0;
278  for (uint16_t iAPV = firstAPV; iAPV < firstAPV + digis.size() / nTotStripsPerAPV; ++iAPV) {
279  if (!badAPVs_[iAPV]) {
280  float MeanAPVCM = meanCM_;
281  if (useRealMeanCM_ && (std::end(meanCMmap_) != itCMMap))
282  MeanAPVCM = itCMMap->second[iAPV];
283 
284  for (uint16_t istrip = iAPV * nTotStripsPerAPV; istrip < (iAPV + 1) * nTotStripsPerAPV; ++istrip) {
285  const auto fs = static_cast<int>(digis[istrip - firstAPV * nTotStripsPerAPV]);
286  if (!qualityHandle->IsStripBad(detQualityRange, istrip)) {
287  qualityCount++;
288  if (std::abs(fs - MeanAPVCM) > static_cast<int>(deviation_)) {
289  devCount++;
290  minstrip = std::min(fs, minstrip);
291  }
292  }
293  }
294  if (devCount > fraction_ * qualityCount) {
295  apvFlags_[iAPV] = restoreAlgo_; //specify any algo to make the restore
296  nAPVflagged++;
297  }
298  }
299  }
300  return nAPVflagged;
301 }
302 
303 inline uint16_t SiStripAPVRestorer::nullInspect(uint16_t firstAPV, const digivector_t& digis) {
305 
306  uint16_t nAPVflagged = 0;
307  for (uint16_t iAPV = firstAPV; iAPV < firstAPV + digis.size() / nTotStripsPerAPV; ++iAPV) {
308  if (!badAPVs_[iAPV]) {
309  int zeroCount = 0, qualityCount = 0;
310  for (uint16_t istrip = iAPV * nTotStripsPerAPV; istrip < (iAPV + 1) * nTotStripsPerAPV; ++istrip) {
311  const auto fs = static_cast<int>(digis[istrip - firstAPV * nTotStripsPerAPV]);
312  if (!qualityHandle->IsStripBad(detQualityRange, istrip)) {
313  qualityCount++;
314  if (fs < 1)
315  zeroCount++;
316  }
317  }
318  if (zeroCount > restoreThreshold_ * qualityCount) {
319  apvFlags_[iAPV] = restoreAlgo_; //specify any algo to make the restore
320  nAPVflagged++;
321  }
322  }
323  }
324  return nAPVflagged;
325 }
326 
327 inline uint16_t SiStripAPVRestorer::hybridEmulationInspect(uint16_t firstAPV, const digivector_t& digis) {
328  uint16_t nAPVflagged = 0;
329 
330  auto itCMMap = std::end(meanCMmap_);
331  if (useRealMeanCM_)
332  itCMMap = meanCMmap_.find(detId_);
333 
334  for (uint16_t iAPV = firstAPV; iAPV < firstAPV + digis.size() / nTotStripsPerAPV; ++iAPV) {
335  if (!badAPVs_[iAPV]) {
336  float MeanAPVCM = meanCM_;
337  if (useRealMeanCM_ && (std::end(meanCMmap_) != itCMMap))
338  MeanAPVCM = itCMMap->second[iAPV];
339 
340  const float DeltaCM = median_[iAPV] - (MeanAPVCM + 1024) / 2;
341  if ((DeltaCM < 0) && (std::abs(DeltaCM) > deltaCMThreshold_ / 2)) {
342  apvFlags_[iAPV] = "HybridEmulation";
343  ++nAPVflagged;
344  }
345  }
346  }
347  return nAPVflagged;
348 }
349 
350 // Restore method implementations
351 
353  uint16_t firstAPV,
354  float median,
355  digivector_t& digis) {
356  digivector_t baseline(size_t(nTotStripsPerAPV), 0);
357 
358  //============================= Find Flat Regions & Interpolating the baseline & subtracting the baseline =================
359  if (!smoothedMaps_.empty()) {
360  auto itSmootedMap = smoothedMaps_.find(apvN);
361  baselineFollower(itSmootedMap->second, baseline, median);
362  } else {
363  //median=0;
364  digivector_t singleAPVdigi;
365  singleAPVdigi.reserve(nTotStripsPerAPV);
366  std::copy_n(
367  std::begin(digis) + nTotStripsPerAPV * (apvN - firstAPV), nTotStripsPerAPV, std::back_inserter(singleAPVdigi));
368  digimap_t smoothedpoints;
369  flatRegionsFinder(singleAPVdigi, smoothedpoints, apvN);
370  baselineFollower(smoothedpoints, baseline, median);
371  }
372 
374  if (checkBaseline(baseline))
375  apvFlagsBoolOverride_[apvN] = true;
376  }
377 
378  //============================= subtracting the baseline =============================================
379  for (int16_t itStrip = 0; itStrip < nTotStripsPerAPV; ++itStrip) {
380  digis[nTotStripsPerAPV * (apvN - firstAPV) + itStrip] -= baseline[itStrip] - median;
381  }
382 
383  //============================= storing baseline to the map =============================================
384  baselineMap_.emplace_hint(baselineMap_.end(), apvN, std::move(baseline));
385 }
386 
387 inline void SiStripAPVRestorer::flatRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis) {
388  digivector_t baseline(size_t(nTotStripsPerAPV), 150);
389  baseline[0] = 0;
390  baseline[nTotStripsPerAPV - 1] = 0;
391 
392  for (int16_t itStrip = 0; itStrip < nTotStripsPerAPV; ++itStrip) {
393  digis[nTotStripsPerAPV * (apvN - firstAPV) + itStrip] = baseline[itStrip];
394  }
395  baselineMap_.emplace_hint(baselineMap_.end(), apvN, std::move(baseline));
396 }
397 
398 // Baseline calculation implementation =========================================
399 bool inline SiStripAPVRestorer::flatRegionsFinder(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN) {
400  smoothedpoints.clear();
401 
403 
404  //============================= Height above local minimum ===============================
405  std::vector<float> adcsLocalMinSubtracted(size_t(nTotStripsPerAPV), 0);
406  for (uint32_t istrip = 0; istrip < nTotStripsPerAPV; ++istrip) {
407  float localmin = 999.9;
408  for (uint16_t jstrip = std::max(0, static_cast<int>(istrip - nSmooth_ / 2));
409  jstrip < std::min((int)nTotStripsPerAPV, static_cast<int>(istrip + nSmooth_ / 2));
410  ++jstrip) {
411  const float nextvalue = adcs[jstrip];
412  if (nextvalue < localmin)
413  localmin = nextvalue;
414  }
415  adcsLocalMinSubtracted[istrip] = adcs[istrip] - localmin;
416  }
417 
418  //============================= Find regions with stable slopes ========================
419  digimap_t consecpoints;
420  std::vector<uint16_t> nConsStrip;
421  //Creating maps with all the neighborhood strip and putting in a nCosntStip vector how many we have
422  uint16_t consecStrips = 0;
423  for (uint32_t istrip = 0; istrip < nTotStripsPerAPV; ++istrip) {
424  const int16_t adc = adcs[istrip];
425  //if( adcsLocalMinSubtracted[istrip] < nSigmaNoiseDerTh_ * (float)noiseHandle->getNoise(istrip+apvN*128,detNoiseRange) && (adc - median) < hitStripThreshold_){
426  if (adcsLocalMinSubtracted[istrip] <
427  nSigmaNoiseDerTh_ * (float)noiseHandle->getNoiseFast(istrip + apvN * nTotStripsPerAPV, detNoiseRange)) {
428  consecpoints.emplace_hint(consecpoints.end(), istrip, adc);
429  ++consecStrips;
430  } else if (consecStrips > 0) {
431  nConsStrip.push_back(consecStrips);
432  consecStrips = 0;
433  }
434  }
435  //to cope with the last flat region of the APV
436  if (consecStrips > 0)
437  nConsStrip.push_back(consecStrips);
438 
439  //removing from the map the fist and last points in wide flat regions and erasing from the map too small regions
440  auto itConsecpoints = consecpoints.begin();
441  float MinSmoothValue = 20000., MaxSmoothValue = 0.;
442  for (const auto consecStrips : nConsStrip) {
443  if (consecStrips >= consecThreshold_) {
444  ++itConsecpoints; //skipping first point
445  uint16_t nFirstStrip = itConsecpoints->first;
446  uint16_t nLastStrip;
447  float smoothValue = 0.0;
448  float stripCount = 1;
449  for (uint16_t n = 0; n < consecStrips - 2; ++n) {
450  smoothValue += itConsecpoints->second;
451  if (stripCount == consecThreshold_) {
452  smoothValue /= (float)stripCount;
453  nLastStrip = nFirstStrip + stripCount - 1;
454  smoothedpoints.emplace_hint(smoothedpoints.end(), nFirstStrip, smoothValue);
455  smoothedpoints.emplace_hint(smoothedpoints.end(), nLastStrip, smoothValue);
456  if (smoothValue > MaxSmoothValue)
457  MaxSmoothValue = smoothValue;
458  if (smoothValue < MinSmoothValue)
459  MinSmoothValue = smoothValue;
460  nFirstStrip = nLastStrip + 1;
461  smoothValue = 0;
462  stripCount = 0;
463  }
464  ++stripCount;
465  ++itConsecpoints;
466  }
467  ++itConsecpoints; //and putting the pointer to the new seies of point
468 
469  if (stripCount > 1) {
470  //if(smoothValue>0){
471  --stripCount;
472  smoothValue /= static_cast<float>(stripCount);
473  nLastStrip = nFirstStrip + stripCount - 1;
474  smoothedpoints.emplace_hint(smoothedpoints.end(), nFirstStrip, smoothValue);
475  smoothedpoints.emplace_hint(smoothedpoints.end(), nLastStrip, smoothValue);
476  if (smoothValue > MaxSmoothValue)
477  MaxSmoothValue = smoothValue;
478  if (smoothValue < MinSmoothValue)
479  MinSmoothValue = smoothValue;
480  }
481  } else {
482  for (int n = 0; n < consecStrips; ++n)
483  ++itConsecpoints;
484  }
485  }
486 
487  if ((MaxSmoothValue - MinSmoothValue) > distortionThreshold_) {
489  baselineCleaner(adcs, smoothedpoints, apvN);
490  return false;
491  }
492  return true;
493 }
494 
495 void inline SiStripAPVRestorer::baselineCleaner(const digivector_t& adcs, digimap_t& smoothedpoints, uint16_t apvN) {
496  if (cleaningSequence_ == 0) { //default sequence used up to now
497  cleaner_HighSlopeChecker(smoothedpoints);
498  cleaner_LocalMinimumAdder(adcs, smoothedpoints, apvN);
499  } else if (cleaningSequence_ == 1) {
500  cleaner_LocalMinimumAdder(adcs, smoothedpoints, apvN);
501  cleaner_HighSlopeChecker(smoothedpoints);
502  cleaner_MonotonyChecker(smoothedpoints);
503  } else if (cleaningSequence_ == 2) {
504  cleaner_HighSlopeChecker(smoothedpoints);
505  } else if (cleaningSequence_ == 3) {
506  cleaner_LocalMinimumAdder(adcs, smoothedpoints, apvN);
507  cleaner_HighSlopeChecker(smoothedpoints);
508  } else {
509  cleaner_HighSlopeChecker(smoothedpoints);
510  cleaner_LocalMinimumAdder(adcs, smoothedpoints, apvN);
511  }
512 }
513 
515  //Removing points without monotony
516  //--------------------------------------------------------------------------------------------------
517  if (smoothedpoints.size() < 3)
518  return;
519 
520  auto it = smoothedpoints.begin();
521  while (smoothedpoints.size() > 3 && it != --(--smoothedpoints.end())) { //while we are not at the last point
522  // get info about current and next points
523  auto it2 = it;
524  const float adc1 = it2->second;
525  const float adc2 = (++it2)->second;
526  const float adc3 = (++it2)->second;
527 
528  if ((adc2 - adc1) > hitStripThreshold_ && (adc2 - adc3) > hitStripThreshold_) {
529  smoothedpoints.erase(--it2);
530  } else {
531  ++it;
532  }
533  }
534 }
535 
537  digimap_t& smoothedpoints,
538  uint16_t apvN) {
540  //inserting extra point is case of local minimum
541  //--------------------------------------------------------------------------------------------------
542  // these should be reset now for the point-insertion that follows
543 
544  if (smoothedpoints.size() >= 2) {
545  for (auto it = smoothedpoints.begin(); it != --smoothedpoints.end(); ++it) {
546  auto itNext = it;
547  ++itNext;
548  const float strip1 = it->first;
549  const float strip2 = itNext->first;
550  const float adc1 = it->second;
551  const float adc2 = itNext->second;
552  const float m = (adc2 - adc1) / (strip2 - strip1);
553 
554  //2,4
555  if ((strip2 - strip1) > slopeX_ && std::abs(adc1 - adc2) > slopeY_) {
556  float itStrip = 1;
557  float strip = itStrip + strip1;
558  while (strip < strip2) {
559  const float adc = adcs[strip];
560  const auto noise =
561  static_cast<float>(noiseHandle->getNoiseFast(strip + apvN * nTotStripsPerAPV, detNoiseRange));
562  if (adc < (adc1 + m * itStrip - 2 * noise)) {
563  smoothedpoints.emplace_hint(itNext, strip, adc);
564  ++it;
565  }
566  ++itStrip;
567  ++strip;
568  }
569  }
570  }
571  }
572 
573  const uint16_t firstStripFlat = smoothedpoints.begin()->first;
574  const uint16_t lastStripFlat = (--smoothedpoints.end())->first;
575  const int16_t firstStripFlatADC = smoothedpoints.begin()->second;
576  const int16_t lastStripFlatADC = (--smoothedpoints.end())->second;
577  if (firstStripFlat > 3) {
578  auto it = smoothedpoints.begin();
579  float strip = 0;
580  while (strip < firstStripFlat) {
581  const float adc = adcs[strip];
582  const auto noise = static_cast<float>(noiseHandle->getNoiseFast(strip + apvN * nTotStripsPerAPV, detNoiseRange));
583  if (adc < (firstStripFlatADC - 2 * noise)) {
584  smoothedpoints.emplace_hint(it, strip, adc);
585  ++it;
586  }
587  ++strip;
588  }
589  }
590  if (lastStripFlat < (nTotStripsPerAPV - 3)) {
591  float strip = lastStripFlat + 1;
592  while (strip < nTotStripsPerAPV) {
593  const float adc = adcs[strip];
594  const auto noise = static_cast<float>(noiseHandle->getNoiseFast(strip + apvN * nTotStripsPerAPV, detNoiseRange));
595  if (adc < (lastStripFlatADC - 2 * noise)) {
596  smoothedpoints.emplace_hint(smoothedpoints.end(), strip, adc);
597  }
598  ++strip;
599  }
600  }
601 }
602 
604  //Removing points in the slope is too high
605  //----------------------------------------------------------------------------
606  if (smoothedpoints.size() < 4)
607  return;
608 
609  auto it = smoothedpoints.begin();
610  while (smoothedpoints.size() > 2 && it != --smoothedpoints.end()) { //while we are not at the last point
611  //if(smoothedpoints.size() <2) break;
612  // get info about current and next points
613  auto itNext = it;
614  ++itNext;
615  const float strip1 = it->first;
616  const float strip2 = itNext->first;
617  const float adc1 = it->second;
618  const float adc2 = itNext->second;
619  const float m = (adc2 - adc1) / (strip2 - strip1);
620 
621  if (m > 2) { // in case of large positive slope, remove next point and try again from same current point
622  smoothedpoints.erase(itNext);
623  } else if (m < -2) { // in case of large negative slope, remove current point and either...
624  // move to next point if we have reached the beginning (post-increment to avoid invalidating pointer during erase) or...
625  if (it == smoothedpoints.begin())
626  smoothedpoints.erase(it++);
627  // try again from the previous point if we have not reached the beginning
628  else
629  smoothedpoints.erase(it--);
630  } else { // in case of a flat enough slope, continue on to the next point
631  ++it;
632  }
633  }
634 }
635 
636 void inline SiStripAPVRestorer::baselineFollower(const digimap_t& smoothedpoints,
637  digivector_t& baseline,
638  float median) {
639  //if not enough points
640  if (smoothedpoints.size() < minStripsToFit_) {
641  baseline.insert(baseline.begin(), nTotStripsPerAPV, median);
642  } else {
643  baseline.assign(size_t(nTotStripsPerAPV), 0);
644 
645  const auto lastIt = --smoothedpoints.end();
646  const uint16_t firstStripFlat = smoothedpoints.begin()->first;
647  const uint16_t lastStripFlat = lastIt->first;
648  const int16_t firstStripFlatADC = smoothedpoints.begin()->second;
649  const int16_t lastStripFlatADC = lastIt->second;
650 
651  //adding here the costant line at the extremities
652  baseline.erase(baseline.begin(), baseline.begin() + firstStripFlat);
653  baseline.insert(baseline.begin(), firstStripFlat, firstStripFlatADC);
654 
655  baseline.erase(baseline.begin() + lastStripFlat, baseline.end());
656  baseline.insert(baseline.end(), nTotStripsPerAPV - lastStripFlat, lastStripFlatADC);
657 
658  //IMPORTANT: the itSmoothedpointsEnd should be at least smaller than smoothedpoints.end() -1
659  for (auto itSmoothedpoints = smoothedpoints.begin(); itSmoothedpoints != (--smoothedpoints.end());
660  ++itSmoothedpoints) {
661  auto itNextSmoothedpoints = itSmoothedpoints;
662  ++itNextSmoothedpoints;
663  const float strip1 = itSmoothedpoints->first;
664  const float strip2 = itNextSmoothedpoints->first;
665  const float adc1 = itSmoothedpoints->second;
666  const float adc2 = itNextSmoothedpoints->second;
667 
668  baseline[strip1] = adc1;
669  baseline[strip2] = adc2;
670  const float m = (adc2 - adc1) / (strip2 - strip1);
671  uint16_t itStrip = strip1 + 1;
672  float stripadc = adc1 + m;
673  while (itStrip < strip2) {
674  baseline[itStrip] = stripadc;
675  ++itStrip;
676  stripadc += m;
677  }
678  }
679  }
680 }
681 
682 bool SiStripAPVRestorer::checkBaseline(const digivector_t& baseline) const {
683  // The Savitzky-Golay (S-G) filter of any left length nL, right
684  // length nR, and order m, and with an optional opt equals the
685  // derivative order (0 for the magnitude, 1 for the first
686  // derivative, etc.) can be calculated using the following
687  // Mathematica code:
688  //
689  // SavitzkyGolay[m_?IntegerQ, {nL_?IntegerQ, nR_?IntegerQ},
690  // opt___] := Module[
691  // {a, d},
692  // d = If[opt === Null, 0, If[IntegerQ[opt], opt, 0]];
693  // a = Table[
694  // If[i == 0 && j == 0, 1, i^j], {i, -nL, nR}, {j, 0,
695  // m}]; (Inverse[Transpose[a].a].Transpose[a])[[d + 1]]];
696  //
697  // The following coefficients can be then calculated by:
698  //
699  // N[Join[Table[SavitzkyGolay[2, {k, 16}], {k, 0, 16}],
700  // Table[SavitzkyGolay[2, {16, k}], {k, 15, 0, -1}]]]
701 
702  // nLR = max(nL, nR)
703  static const size_t savitzky_golay_n_l_r = 16;
704  static const float savitzky_golay_coefficient[2 * savitzky_golay_n_l_r + 1][2 * savitzky_golay_n_l_r + 1] = {
705  {0.422085,
706  0.325077,
707  0.23839,
708  0.162023,
709  0.0959752,
710  0.0402477,
711  -0.00515996,
712  -0.0402477,
713  -0.0650155,
714  -0.0794634,
715  -0.0835913,
716  -0.0773994,
717  -0.0608875,
718  -0.0340557,
719  0.00309598,
720  0.0505676,
721  0.108359},
722  {0.315789,
723  0.254902,
724  0.19969,
725  0.150155,
726  0.106295,
727  0.0681115,
728  0.0356037,
729  0.00877193,
730  -0.0123839,
731  -0.0278638,
732  -0.0376677,
733  -0.0417957,
734  -0.0402477,
735  -0.0330237,
736  -0.0201238,
737  -0.00154799,
738  0.0227038,
739  0.0526316},
740  {0.234586,
741  0.198496,
742  0.165207,
743  0.134719,
744  0.107032,
745  0.0821465,
746  0.0600619,
747  0.0407784,
748  0.024296,
749  0.0106148,
750  -0.000265369,
751  -0.00834439,
752  -0.0136223,
753  -0.0160991,
754  -0.0157747,
755  -0.0126493,
756  -0.00672269,
757  0.00200501,
758  0.0135338},
759  {0.172078, 0.153076, 0.135099, 0.118148, 0.102221, 0.0873206, 0.073445,
760  0.0605947, 0.0487697, 0.0379699, 0.0281955, 0.0194463, 0.0117225, 0.00502392,
761  -0.000649351, -0.00529733, -0.00892003, -0.0115174, -0.0130895, -0.0136364},
762  {0.123659, 0.116431, 0.109144, 0.101798, 0.0943921, 0.0869268, 0.0794021,
763  0.0718179, 0.0641743, 0.0564712, 0.0487087, 0.0408868, 0.0330054, 0.0250646,
764  0.0170644, 0.00900473, 0.000885613, -0.00729294, -0.0155309, -0.0238283, -0.0321852},
765  {0.0859684, 0.0868154, 0.0869565, 0.0863919, 0.0851214, 0.0831451, 0.080463, 0.0770751,
766  0.0729814, 0.0681818, 0.0626765, 0.0564653, 0.0495483, 0.0419255, 0.0335968, 0.0245624,
767  0.0148221, 0.00437606, -0.00677583, -0.0186335, -0.0311971, -0.0444664},
768  {0.0565217, 0.0628458, 0.0680971, 0.0722756, 0.0753811, 0.0774139, 0.0783738, 0.0782609,
769  0.0770751, 0.0748165, 0.071485, 0.0670807, 0.0616036, 0.0550536, 0.0474308, 0.0387352,
770  0.0289667, 0.0181254, 0.00621118, -0.00677583, -0.0208357, -0.0359684, -0.0521739},
771  {0.0334615, 0.0434281, 0.0521329, 0.0595759, 0.0657571, 0.0706765, 0.0743341, 0.07673,
772  0.0778641, 0.0777364, 0.0763469, 0.0736957, 0.0697826, 0.0646078, 0.0581712, 0.0504728,
773  0.0415126, 0.0312907, 0.0198069, 0.00706142, -0.00694588, -0.022215, -0.0387458, -0.0565385},
774  {0.0153846, 0.0276923, 0.0386622, 0.0482943, 0.0565886, 0.0635452, 0.0691639, 0.0734448, 0.076388,
775  0.0779933, 0.0782609, 0.0771906, 0.0747826, 0.0710368, 0.0659532, 0.0595318, 0.0517726, 0.0426756,
776  0.0322408, 0.0204682, 0.00735786, -0.0070903, -0.0228763, -0.04, -0.0584615},
777  {0.001221, 0.0149451, 0.027326, 0.0383639, 0.0480586, 0.0564103, 0.0634188, 0.0690842, 0.0734066,
778  0.0763858, 0.078022, 0.078315, 0.077265, 0.0748718, 0.0711355, 0.0660562, 0.0596337, 0.0518681,
779  0.0427595, 0.0323077, 0.0205128, 0.00737485, -0.00710623, -0.0229304, -0.0400977, -0.0586081},
780  {-0.00985222, 0.00463138, 0.0178098, 0.029683, 0.0402509, 0.0495137, 0.0574713, 0.0641236, 0.0694708,
781  0.0735127, 0.0762494, 0.0776809, 0.0778073, 0.0766284, 0.0741442, 0.0703549, 0.0652604, 0.0588607,
782  0.0511557, 0.0421456, 0.0318302, 0.0202097, 0.0072839, -0.00694708, -0.0224833, -0.0393247, -0.0574713},
783  {-0.0184729, -0.00369458, 0.00984169, 0.0221359, 0.0331881, 0.0429982, 0.0515662,
784  0.0588923, 0.0649762, 0.0698181, 0.073418, 0.0757758, 0.0768915, 0.0767652,
785  0.0753968, 0.0727864, 0.0689339, 0.0638394, 0.0575028, 0.0499242, 0.0411035,
786  0.0310408, 0.019736, 0.00718917, -0.00659972, -0.0216307, -0.0379037, -0.0554187},
787  {-0.025139, -0.0103925, 0.00318873, 0.0156046, 0.0268552, 0.0369405, 0.0458605, 0.0536151,
788  0.0602045, 0.0656285, 0.0698872, 0.0729806, 0.0749086, 0.0756714, 0.0752688, 0.0737009,
789  0.0709677, 0.0670692, 0.0620054, 0.0557763, 0.0483818, 0.039822, 0.0300969, 0.0192065,
790  0.0071508, -0.00607024, -0.0204566, -0.0360083, -0.0527253},
791  {-0.0302419, -0.0157536, -0.00234785, 0.00997537, 0.021216, 0.0313741, 0.0404497, 0.0484427,
792  0.0553532, 0.0611811, 0.0659264, 0.0695892, 0.0721695, 0.0736672, 0.0740823, 0.0734149,
793  0.0716649, 0.0688324, 0.0649174, 0.0599198, 0.0538396, 0.0466769, 0.0384316, 0.0291038,
794  0.0186934, 0.00720046, -0.00537502, -0.0190331, -0.0337736, -0.0495968},
795  {-0.0340909, -0.0200147, -0.006937, 0.00514208, 0.0162226, 0.0263045, 0.0353878, 0.0434725,
796  0.0505587, 0.0566463, 0.0617353, 0.0658257, 0.0689175, 0.0710107, 0.0721054, 0.0722014,
797  0.0712989, 0.0693978, 0.0664981, 0.0625999, 0.057703, 0.0518076, 0.0449135, 0.0370209,
798  0.0281297, 0.01824, 0.0073516, -0.00453534, -0.0174209, -0.031305, -0.0461877},
799  {-0.0369318, -0.0233688, -0.0107221, 0.00100806, 0.0118218, 0.0217192, 0.0307001, 0.0387647,
800  0.0459128, 0.0521444, 0.0574597, 0.0618585, 0.0653409, 0.0679069, 0.0695565, 0.0702896,
801  0.0701063, 0.0690066, 0.0669905, 0.0640579, 0.0602089, 0.0554435, 0.0497617, 0.0431635,
802  0.0356488, 0.0272177, 0.0178702, 0.0076063, -0.00357405, -0.0156708, -0.028684, -0.0426136},
803  {-0.038961, -0.025974, -0.0138249, -0.00251362, 0.00795978, 0.0175953, 0.026393, 0.0343527, 0.0414747,
804  0.0477587, 0.0532049, 0.0578132, 0.0615836, 0.0645161, 0.0666108, 0.0678676, 0.0682866, 0.0678676,
805  0.0666108, 0.0645161, 0.0615836, 0.0578132, 0.0532049, 0.0477587, 0.0414747, 0.0343527, 0.026393,
806  0.0175953, 0.00795978, -0.00251362, -0.0138249, -0.025974, -0.038961},
807  {-0.0426136, -0.028684, -0.0156708, -0.00357405, 0.0076063, 0.0178702, 0.0272177, 0.0356488,
808  0.0431635, 0.0497617, 0.0554435, 0.0602089, 0.0640579, 0.0669905, 0.0690066, 0.0701063,
809  0.0702896, 0.0695565, 0.0679069, 0.0653409, 0.0618585, 0.0574597, 0.0521444, 0.0459128,
810  0.0387647, 0.0307001, 0.0217192, 0.0118218, 0.00100806, -0.0107221, -0.0233688, -0.0369318},
811  {-0.0461877, -0.031305, -0.0174209, -0.00453534, 0.0073516, 0.01824, 0.0281297, 0.0370209,
812  0.0449135, 0.0518076, 0.057703, 0.0625999, 0.0664981, 0.0693978, 0.0712989, 0.0722014,
813  0.0721054, 0.0710107, 0.0689175, 0.0658257, 0.0617353, 0.0566463, 0.0505587, 0.0434725,
814  0.0353878, 0.0263045, 0.0162226, 0.00514208, -0.006937, -0.0200147, -0.0340909},
815  {-0.0495968, -0.0337736, -0.0190331, -0.00537502, 0.00720046, 0.0186934, 0.0291038, 0.0384316,
816  0.0466769, 0.0538396, 0.0599198, 0.0649174, 0.0688324, 0.0716649, 0.0734149, 0.0740823,
817  0.0736672, 0.0721695, 0.0695892, 0.0659264, 0.0611811, 0.0553532, 0.0484427, 0.0404497,
818  0.0313741, 0.021216, 0.00997537, -0.00234785, -0.0157536, -0.0302419},
819  {-0.0527253, -0.0360083, -0.0204566, -0.00607024, 0.0071508, 0.0192065, 0.0300969, 0.039822,
820  0.0483818, 0.0557763, 0.0620054, 0.0670692, 0.0709677, 0.0737009, 0.0752688, 0.0756714,
821  0.0749086, 0.0729806, 0.0698872, 0.0656285, 0.0602045, 0.0536151, 0.0458605, 0.0369405,
822  0.0268552, 0.0156046, 0.00318873, -0.0103925, -0.025139},
823  {-0.0554187, -0.0379037, -0.0216307, -0.00659972, 0.00718917, 0.019736, 0.0310408,
824  0.0411035, 0.0499242, 0.0575028, 0.0638394, 0.0689339, 0.0727864, 0.0753968,
825  0.0767652, 0.0768915, 0.0757758, 0.073418, 0.0698181, 0.0649762, 0.0588923,
826  0.0515662, 0.0429982, 0.0331881, 0.0221359, 0.00984169, -0.00369458, -0.0184729},
827  {-0.0574713, -0.0393247, -0.0224833, -0.00694708, 0.0072839, 0.0202097, 0.0318302, 0.0421456, 0.0511557,
828  0.0588607, 0.0652604, 0.0703549, 0.0741442, 0.0766284, 0.0778073, 0.0776809, 0.0762494, 0.0735127,
829  0.0694708, 0.0641236, 0.0574713, 0.0495137, 0.0402509, 0.029683, 0.0178098, 0.00463138, -0.00985222},
830  {-0.0586081, -0.0400977, -0.0229304, -0.00710623, 0.00737485, 0.0205128, 0.0323077, 0.0427595, 0.0518681,
831  0.0596337, 0.0660562, 0.0711355, 0.0748718, 0.077265, 0.078315, 0.078022, 0.0763858, 0.0734066,
832  0.0690842, 0.0634188, 0.0564103, 0.0480586, 0.0383639, 0.027326, 0.0149451, 0.001221},
833  {-0.0584615, -0.04, -0.0228763, -0.0070903, 0.00735786, 0.0204682, 0.0322408, 0.0426756, 0.0517726,
834  0.0595318, 0.0659532, 0.0710368, 0.0747826, 0.0771906, 0.0782609, 0.0779933, 0.076388, 0.0734448,
835  0.0691639, 0.0635452, 0.0565886, 0.0482943, 0.0386622, 0.0276923, 0.0153846},
836  {-0.0565385, -0.0387458, -0.022215, -0.00694588, 0.00706142, 0.0198069, 0.0312907, 0.0415126,
837  0.0504728, 0.0581712, 0.0646078, 0.0697826, 0.0736957, 0.0763469, 0.0777364, 0.0778641,
838  0.07673, 0.0743341, 0.0706765, 0.0657571, 0.0595759, 0.0521329, 0.0434281, 0.0334615},
839  {-0.0521739, -0.0359684, -0.0208357, -0.00677583, 0.00621118, 0.0181254, 0.0289667, 0.0387352,
840  0.0474308, 0.0550536, 0.0616036, 0.0670807, 0.071485, 0.0748165, 0.0770751, 0.0782609,
841  0.0783738, 0.0774139, 0.0753811, 0.0722756, 0.0680971, 0.0628458, 0.0565217},
842  {-0.0444664, -0.0311971, -0.0186335, -0.00677583, 0.00437606, 0.0148221, 0.0245624, 0.0335968,
843  0.0419255, 0.0495483, 0.0564653, 0.0626765, 0.0681818, 0.0729814, 0.0770751, 0.080463,
844  0.0831451, 0.0851214, 0.0863919, 0.0869565, 0.0868154, 0.0859684},
845  {-0.0321852, -0.0238283, -0.0155309, -0.00729294, 0.000885613, 0.00900473, 0.0170644,
846  0.0250646, 0.0330054, 0.0408868, 0.0487087, 0.0564712, 0.0641743, 0.0718179,
847  0.0794021, 0.0869268, 0.0943921, 0.101798, 0.109144, 0.116431, 0.123659},
848  {-0.0136364, -0.0130895, -0.0115174, -0.00892003, -0.00529733, -0.000649351, 0.00502392,
849  0.0117225, 0.0194463, 0.0281955, 0.0379699, 0.0487697, 0.0605947, 0.073445,
850  0.0873206, 0.102221, 0.118148, 0.135099, 0.153076, 0.172078},
851  {0.0135338,
852  0.00200501,
853  -0.00672269,
854  -0.0126493,
855  -0.0157747,
856  -0.0160991,
857  -0.0136223,
858  -0.00834439,
859  -0.000265369,
860  0.0106148,
861  0.024296,
862  0.0407784,
863  0.0600619,
864  0.0821465,
865  0.107032,
866  0.134719,
867  0.165207,
868  0.198496,
869  0.234586},
870  {0.0526316,
871  0.0227038,
872  -0.00154799,
873  -0.0201238,
874  -0.0330237,
875  -0.0402477,
876  -0.0417957,
877  -0.0376677,
878  -0.0278638,
879  -0.0123839,
880  0.00877193,
881  0.0356037,
882  0.0681115,
883  0.106295,
884  0.150155,
885  0.19969,
886  0.254902,
887  0.315789},
888  {0.108359,
889  0.0505676,
890  0.00309598,
891  -0.0340557,
892  -0.0608875,
893  -0.0773994,
894  -0.0835913,
895  -0.0794634,
896  -0.0650155,
897  -0.0402477,
898  -0.00515996,
899  0.0402477,
900  0.0959752,
901  0.162023,
902  0.23839,
903  0.325077,
904  0.422085}};
905 
906  float filtered_baseline[nTotStripsPerAPV];
907  float filtered_baseline_derivative[nTotStripsPerAPV - 1];
908 
909  // Zero filtered_baseline
910  memset(filtered_baseline, 0, nTotStripsPerAPV * sizeof(float));
911  // Filter the left edge using (nL, nR) = (0, 16) .. (15, 16) S-G
912  // filters
913  for (size_t i = 0; i < savitzky_golay_n_l_r; i++) {
914  for (size_t j = 0; j < savitzky_golay_n_l_r + 1 + i; j++) {
915  filtered_baseline[i] += savitzky_golay_coefficient[i][j] * baseline[j];
916  }
917  }
918  // Filter the middle section using the (nL, nR) = (16, 16) S-G
919  // filter, while taking advantage of the symmetry to save 16
920  // multiplications.
921  for (size_t i = savitzky_golay_n_l_r; i < nTotStripsPerAPV - savitzky_golay_n_l_r; i++) {
922  filtered_baseline[i] = savitzky_golay_coefficient[savitzky_golay_n_l_r][savitzky_golay_n_l_r] * baseline[i];
923  for (size_t j = 0; j < savitzky_golay_n_l_r; j++) {
924  filtered_baseline[i] += savitzky_golay_coefficient[savitzky_golay_n_l_r][j] *
925  (baseline[i + j - savitzky_golay_n_l_r] + baseline[i - j + savitzky_golay_n_l_r]);
926  }
927 #if 0
928  // Test that the indexing above is correct
929  float test = 0;
930  for (size_t j = 0; j < 2 * savitzky_golay_n_l_r + 1; j++) {
931  test +=
932  savitzky_golay_coefficient[savitzky_golay_n_l_r][j] *
933  baseline[i + j - savitzky_golay_n_l_r];
934  }
935  // test == filtered_baseline[i] should hold now
936 #endif
937  }
938  // Filter the right edge using (nL, nR) = (16, 15) .. (16, 0) S-G
939  // filters
940  for (size_t i = nTotStripsPerAPV - savitzky_golay_n_l_r; i < nTotStripsPerAPV; i++) {
941  for (size_t j = 0; j < nTotStripsPerAPV - i + savitzky_golay_n_l_r; j++) {
942  filtered_baseline[i] += savitzky_golay_coefficient[2 * savitzky_golay_n_l_r + i + 1 - nTotStripsPerAPV][j] *
943  baseline[i + j - savitzky_golay_n_l_r];
944  }
945  }
946  // In lieu of a spearate S-G derivative filter, the finite
947  // difference is used here (since the output is sufficiently
948  // smooth).
949  for (size_t i = 0; i < (nTotStripsPerAPV - 1); i++) {
950  filtered_baseline_derivative[i] = filtered_baseline[i + 1] - filtered_baseline[i];
951  }
952 
953  // Calculate the maximum deviation between filtered and unfiltered
954  // baseline, plus the sum square of the derivative.
955 
956  double filtered_baseline_max = 0;
957  double filtered_baseline_derivative_sum_square = 0;
958 
959  for (size_t i = 0; i < nTotStripsPerAPV; i++) {
960  const double d = filtered_baseline[i] - baseline[i];
961 
962  filtered_baseline_max = std::max(filtered_baseline_max, static_cast<double>(fabs(d)));
963  }
964  for (size_t i = 0; i < (nTotStripsPerAPV - 1); i++) {
965  filtered_baseline_derivative_sum_square += filtered_baseline_derivative[i] * filtered_baseline_derivative[i];
966  }
967 
968 #if 0
969  std::cerr << __FILE__ << ':' << __LINE__ << ": "
970  << filtered_baseline_max << ' '
971  << filtered_baseline_derivative_sum_square << std::endl;
972 #endif
973 
974  // Apply the cut
975  return !(filtered_baseline_max >= filteredBaselineMax_ ||
976  filtered_baseline_derivative_sum_square >= filteredBaselineDerivativeSumSquare_);
977 }
978 
979 //Other methods implementation ==============================================
980 //==========================================================================
981 
983  if (useRealMeanCM_) {
985  iEvent.getByLabel("siStripDigis", "VirginRaw", input);
987  } else {
989  iEvent.getByLabel("MEANAPVCM", inputCM);
990  createCMMapCMstored(*inputCM);
991  }
992 }
994  meanCMmap_.clear();
995  for (const auto& rawDigis : input) {
996  const auto detPedestalRange = pedestalHandle->getRange(rawDigis.id);
997  const size_t nAPVs = rawDigis.size() / nTotStripsPerAPV;
998  std::vector<float> meanCMDetSet;
999  meanCMDetSet.reserve(nAPVs);
1000  for (uint16_t iAPV = 0; iAPV < nAPVs; ++iAPV) {
1001  uint16_t minPed = nTotStripsPerAPV;
1002  for (uint16_t strip = nTotStripsPerAPV * iAPV; strip < nTotStripsPerAPV * (iAPV + 1); ++strip) {
1003  uint16_t ped = static_cast<uint16_t>(pedestalHandle->getPed(strip, detPedestalRange));
1004  if (ped < minPed)
1005  minPed = ped;
1006  }
1007  meanCMDetSet.push_back(minPed);
1008  }
1009  meanCMmap_.emplace(rawDigis.id, std::move(meanCMDetSet));
1010  }
1011 }
1013  meanCMmap_.clear();
1014  for (const auto& rawDigis : input) {
1015  std::vector<float> meanCMNValue;
1016  meanCMNValue.reserve(rawDigis.size());
1018  std::begin(rawDigis), std::end(rawDigis), std::back_inserter(meanCMNValue), [](SiStripProcessedRawDigi cm) {
1019  return cm.adc();
1020  });
1021  meanCMmap_.emplace(rawDigis.id, std::move(meanCMNValue));
1022  }
1023 }
1024 
1025 //NEW algorithm designed for implementation in FW=============================
1026 //============================================================================
1027 void SiStripAPVRestorer::derivativeFollowerRestore(uint16_t apvN, uint16_t firstAPV, digivector_t& digis) {
1028  digivector_t singleAPVdigi;
1029  singleAPVdigi.clear();
1030  for (int16_t strip = (apvN - firstAPV) * nTotStripsPerAPV; strip < (apvN - firstAPV + 1) * nTotStripsPerAPV; ++strip)
1031  singleAPVdigi.push_back(digis[strip] + 1024);
1032 
1033  digimap_t discontinuities; // it will contain the start and the end of each region in which a greadient is present
1034  discontinuities.clear();
1035 
1036  digimap_t::iterator itdiscontinuities;
1037 
1038  //----Variables of the first part---//
1039 
1040  bool isMinimumAndNoMax = false;
1041  bool isFirstStrip = false;
1042 
1043  //---Variables of the second part---//
1044 
1045  int actualStripADC = 0;
1046  int previousStripADC = 0;
1047 
1048  int greadient = 0;
1049  int maximum_value = 0;
1050  const uint32_t n_high_maximum_cluster = 1025 + 1024;
1051  int high_maximum_cluster = n_high_maximum_cluster;
1052  int number_good_minimum = 0;
1053  int first_gradient = 0;
1054  int strip_first_gradient = 0;
1055  int adc_start_point_cluster_pw = 0;
1056  int auxiliary_end_cluster = 0;
1057  int first_start_cluster_strip = 0;
1058  int first_start_cluster_ADC = 0;
1059  bool isAuxiliary_Minimum = false;
1060  bool isPossible_wrong_minimum = false;
1061  bool isMax = false;
1062 
1063  //----------SECOND PART: CLUSTER FINDING--------//
1064 
1065  for (uint16_t strip = 0; strip < singleAPVdigi.size(); ++strip) {
1066  if (strip == 0) {
1067  actualStripADC = singleAPVdigi[strip];
1068  if (std::abs(singleAPVdigi[strip] - singleAPVdigi[strip + 1]) > gradient_threshold_) {
1069  isFirstStrip = true;
1070  isMinimumAndNoMax = true;
1071  discontinuities.insert(discontinuities.end(), std::pair<int, int>(strip, actualStripADC));
1072  } else if (actualStripADC > (gradient_threshold_ + 1024)) {
1073  discontinuities.insert(discontinuities.end(), std::pair<uint16_t, int16_t>(strip, 0 + 1024));
1074  isMinimumAndNoMax = true;
1075  first_start_cluster_strip = strip;
1076  first_start_cluster_ADC = 1024;
1077  }
1078 
1079  } else {
1080  previousStripADC = actualStripADC;
1081  actualStripADC = singleAPVdigi[strip];
1082  greadient = actualStripADC - previousStripADC;
1083 
1084  if (((greadient > gradient_threshold_) && (!isMax) && (!isMinimumAndNoMax))) {
1085  isMax = false;
1086  if (((std::abs(maximum_value - previousStripADC) < (2 * gradient_threshold_)) &&
1087  discontinuities.size() > 1)) { // add the || with the cluster size
1088  //to verify if the noise do not interfere and is detected fake hits
1089  isPossible_wrong_minimum = true;
1090  itdiscontinuities = --discontinuities.end();
1091  if (discontinuities.size() > 1) {
1092  --itdiscontinuities;
1093  }
1094  strip_first_gradient = itdiscontinuities->first;
1095  adc_start_point_cluster_pw = itdiscontinuities->second;
1096  first_gradient = std::abs(adc_start_point_cluster_pw - singleAPVdigi[strip_first_gradient + 1]);
1097  discontinuities.erase(++itdiscontinuities);
1098  discontinuities.erase(--discontinuities.end());
1099  }
1100 
1101  if ((discontinuities.size() % 2 == 1) && (!discontinuities.empty())) { //&&(no_minimo == 0)
1102  discontinuities.erase(--discontinuities.end());
1103  }
1104 
1105  discontinuities.insert(discontinuities.end(), std::pair<uint16_t, int16_t>(strip - 1, previousStripADC));
1106  isMinimumAndNoMax = true;
1107  maximum_value = 0;
1108  first_start_cluster_strip = strip - 1;
1109  first_start_cluster_ADC = previousStripADC;
1110 
1111  } else if ((!isMax) && ((actualStripADC - previousStripADC < 0) && (isMinimumAndNoMax))) {
1112  isMax = true;
1113  isMinimumAndNoMax = false;
1114  high_maximum_cluster = n_high_maximum_cluster;
1115  if ((previousStripADC > maximum_value) && (discontinuities.size() % 2 == 1))
1116  maximum_value = previousStripADC;
1117  }
1118 
1119  if ((isMax) && (strip < (nTotStripsPerAPV - 2))) {
1120  if (high_maximum_cluster > (std::abs(singleAPVdigi[strip + 1] - actualStripADC))) {
1121  high_maximum_cluster = singleAPVdigi[strip + 1] - actualStripADC;
1122  auxiliary_end_cluster = strip + 2;
1123  } else {
1124  auxiliary_end_cluster = nTotStripsPerAPV - 1;
1125  }
1126  }
1127 
1128  if ((isMax) && ((actualStripADC - previousStripADC) >= 0) && (size_window_ > 0) &&
1129  (strip <= ((nTotStripsPerAPV - 1) - size_window_ - 1))) {
1130  number_good_minimum = 0;
1131  for (uint16_t wintry = 0; wintry <= size_window_; wintry++) {
1132  if (std::abs(singleAPVdigi[strip + wintry] - singleAPVdigi[strip + wintry + 1]) <= last_gradient_)
1133  ++number_good_minimum;
1134  }
1135  --number_good_minimum;
1136 
1137  if (size_window_ != number_good_minimum) {
1138  isMax = false; // not Valid end Point
1139  isMinimumAndNoMax = true;
1140  }
1141 
1142  } else if ((isMax) && (strip > ((nTotStripsPerAPV - 1) - size_window_ - 1))) {
1143  isMax = true;
1144  isAuxiliary_Minimum = true; //for minimums after strip 127-SW-1
1145  }
1146 
1147  if (!discontinuities.empty()) {
1148  itdiscontinuities = --discontinuities.end();
1149  }
1150 
1151  if ((isMax) && (actualStripADC <= first_start_cluster_ADC)) {
1152  if ((std::abs(first_start_cluster_ADC - singleAPVdigi[strip + 2]) > first_gradient) &&
1153  (isPossible_wrong_minimum)) {
1154  discontinuities.erase(itdiscontinuities);
1155  discontinuities.insert(discontinuities.end(),
1156  std::pair<int, int>(strip_first_gradient, adc_start_point_cluster_pw));
1157  discontinuities.insert(discontinuities.end(), std::pair<int, int>(strip, adc_start_point_cluster_pw));
1158  } else {
1159  if ((discontinuities.size() % 2 == 0) && (discontinuities.size() > 1)) {
1160  discontinuities.erase(--discontinuities.end());
1161  }
1162  discontinuities.insert(discontinuities.end(), std::pair<int, int>(strip, actualStripADC));
1163  }
1164  isMax = false;
1165  adc_start_point_cluster_pw = 0;
1166  isPossible_wrong_minimum = false;
1167  strip_first_gradient = 0;
1168  first_start_cluster_strip = 0;
1169  first_start_cluster_ADC = 0;
1170  }
1171 
1172  if ((isMax) && ((actualStripADC - previousStripADC) >= 0) &&
1173  (!isAuxiliary_Minimum)) { //For the end Poit when strip >127-SW-1
1174  if ((std::abs(first_start_cluster_ADC - singleAPVdigi[strip + 1]) > first_gradient) &&
1175  (isPossible_wrong_minimum)) {
1176  discontinuities.erase(itdiscontinuities);
1177  discontinuities.insert(discontinuities.end(),
1178  std::pair<int, int>(strip_first_gradient, adc_start_point_cluster_pw));
1179  }
1180 
1181  if ((discontinuities.size() % 2 == 0) && (discontinuities.size() > 1)) {
1182  discontinuities.erase(--discontinuities.end());
1183  }
1184  discontinuities.insert(discontinuities.end(), std::pair<int, int>(strip - 1, previousStripADC));
1185  isMax = false;
1186  adc_start_point_cluster_pw = 0;
1187  isPossible_wrong_minimum = false;
1188  strip_first_gradient = 0;
1189  first_start_cluster_strip = 0;
1190  first_start_cluster_ADC = 0;
1191  }
1192  }
1193  }
1194 
1195  //----------THIRD PART reconstruction of the event without baseline-------//
1196 
1197  if (!discontinuities.empty()) {
1198  if ((first_start_cluster_strip) == (nTotStripsPerAPV - 1) - 1)
1199  discontinuities.insert(discontinuities.end(),
1200  std::pair<int, int>((nTotStripsPerAPV - 1), first_start_cluster_ADC));
1201  if ((isMax) && (isAuxiliary_Minimum))
1202  discontinuities.insert(discontinuities.end(),
1203  std::pair<int, int>(auxiliary_end_cluster, first_start_cluster_ADC));
1204  }
1205 
1206  if (isFirstStrip) {
1207  itdiscontinuities = discontinuities.begin();
1208  ++itdiscontinuities;
1209  ++itdiscontinuities;
1210  int firstADC = itdiscontinuities->second;
1211  --itdiscontinuities;
1212  itdiscontinuities->second = firstADC;
1213  --itdiscontinuities;
1214  itdiscontinuities->second = firstADC;
1215  }
1216 
1217  if (!discontinuities.empty()) {
1218  itdiscontinuities = discontinuities.begin();
1219  uint16_t firstStrip = itdiscontinuities->first;
1220  int16_t firstADC = itdiscontinuities->second;
1221  ++itdiscontinuities;
1222  uint16_t lastStrip = itdiscontinuities->first;
1223  int16_t secondADC = itdiscontinuities->second;
1224  ++itdiscontinuities;
1225 
1226  for (uint16_t strip = 0; strip < nTotStripsPerAPV; ++strip) {
1227  if (strip > lastStrip && itdiscontinuities != discontinuities.end()) {
1228  firstStrip = itdiscontinuities->first;
1229  firstADC = itdiscontinuities->second;
1230  ++itdiscontinuities;
1231  lastStrip = itdiscontinuities->first;
1232  secondADC = itdiscontinuities->second;
1233  ++itdiscontinuities;
1234  }
1235 
1236  if ((firstStrip <= strip) && (strip <= lastStrip) &&
1237  (0 < (singleAPVdigi[strip] - firstADC -
1238  ((secondADC - firstADC) / (lastStrip - firstStrip)) * (strip - firstStrip) - gradient_threshold_))) {
1239  digis[(apvN - firstAPV) * nTotStripsPerAPV + strip] =
1240  singleAPVdigi[strip] - firstADC -
1241  (((secondADC - firstADC) / (lastStrip - firstStrip)) * (strip - firstStrip));
1242  edm::LogWarning("NoBaseline") << "No baseline " << digis[(apvN - firstAPV) * nTotStripsPerAPV + strip] << "\n";
1243  } else {
1244  digis[(apvN - firstAPV) * nTotStripsPerAPV + strip] = 0;
1245  }
1246  }
1247  }
1248 }
edm::DetSetVector< SiStripRawDigi >
ecalMGPA::adc
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
Definition: EcalMGPASample.h:11
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
input
static const std::string input
Definition: EdmProvDump.cc:48
SiStripAPVRestorer::pedestal_cache_id
uint32_t pedestal_cache_id
Definition: SiStripAPVRestorer.h:89
MessageLogger.h
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
SiStripNoises::getNoiseFast
static float getNoiseFast(const uint16_t &strip, const Range &range)
Definition: SiStripNoises.h:66
SiStripBadStrip::getRange
const Range getRange(const uint32_t detID) const
Definition: SiStripBadStrip.cc:27
SiStripAPVRestorer::deviation_
uint32_t deviation_
Definition: SiStripAPVRestorer.h:113
SiStripAPVRestorer::smoothedMaps_
std::map< uint16_t, digimap_t > smoothedMaps_
Definition: SiStripAPVRestorer.h:100
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
SiStripAPVRestorer::inspectAlgo_
std::string inspectAlgo_
Definition: SiStripAPVRestorer.h:107
SiStripAPVRestorer::filteredBaselineDerivativeSumSquare_
double filteredBaselineDerivativeSumSquare_
Definition: SiStripAPVRestorer.h:128
min
T min(T a, T b)
Definition: MathUtil.h:58
SiStripAPVRestorer::gradient_threshold_
int gradient_threshold_
Definition: SiStripAPVRestorer.h:129
digitizers_cfi.strip
strip
Definition: digitizers_cfi.py:19
SiStripAPVRestorer::pedestalHandle
edm::ESHandle< SiStripPedestals > pedestalHandle
Definition: SiStripAPVRestorer.h:88
SiStripAPVRestorer::loadMeanCMMap
void loadMeanCMMap(const edm::Event &)
Definition: SiStripAPVRestorer.cc:982
SiStripAPVRestorer::nSmooth_
uint32_t nSmooth_
Definition: SiStripAPVRestorer.h:118
SiStripAPVRestorer::baselineMap_
baselinemap_t baselineMap_
Definition: SiStripAPVRestorer.h:101
SiStripAPVRestorer::forceNoRestore_
bool forceNoRestore_
Definition: SiStripAPVRestorer.h:106
SiStripAPVRestorer::baselineFollower
void baselineFollower(const digimap_t &, digivector_t &baseline, float median)
Definition: SiStripAPVRestorer.cc:636
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:215
SiStripProcessedRawDigi::adc
const float & adc() const
Definition: SiStripProcessedRawDigi.h:36
SiStripNoises::getRange
const Range getRange(const uint32_t detID) const
Definition: SiStripNoises.cc:34
SiStripAPVRestorer::nSaturatedStrip_
uint32_t nSaturatedStrip_
Definition: SiStripAPVRestorer.h:115
edm::Handle
Definition: AssociativeIterator.h:50
dqmdumpme.first
first
Definition: dqmdumpme.py:55
SiStripAPVRestorer::init
void init(const edm::EventSetup &es)
Definition: SiStripAPVRestorer.cc:54
SiStripPedestalsRcd
Definition: SiStripCondDataRecords.h:44
SiStripNoisesRcd
Definition: SiStripCondDataRecords.h:40
SiStripAPVRestorer::cleaningSequence_
uint32_t cleaningSequence_
Definition: SiStripAPVRestorer.h:121
SiStripAPVRestorer::minStripsToFit_
uint32_t minStripsToFit_
Definition: SiStripAPVRestorer.h:125
SiStripAPVRestorer::checkBaseline
bool checkBaseline(const std::vector< int16_t > &baseline) const
Definition: SiStripAPVRestorer.cc:682
end
#define end
Definition: vmac.h:39
SiStripAPVRestorer::restoreAlgo_
std::string restoreAlgo_
Definition: SiStripAPVRestorer.h:108
SiStripAPVRestorer::applyBaselineRejection_
bool applyBaselineRejection_
Definition: SiStripAPVRestorer.h:126
SiStripAPVRestorer::cleaner_LocalMinimumAdder
void cleaner_LocalMinimumAdder(const digivector_t &adcs, digimap_t &smoothedpoints, uint16_t apvN)
Definition: SiStripAPVRestorer.cc:536
SiStripAPVRestorer::deltaCMThreshold_
uint32_t deltaCMThreshold_
Definition: SiStripAPVRestorer.h:111
SiStripAPVRestorer::cleaner_HighSlopeChecker
void cleaner_HighSlopeChecker(digimap_t &smoothedpoints)
Definition: SiStripAPVRestorer.cc:603
SiStripBadStrip::Range
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripBadStrip.h:53
test
Definition: SmallWORMDict.h:13
SiStripAPVRestorer::slopeY_
int32_t slopeY_
Definition: SiStripAPVRestorer.h:123
SiStripAPVRestorer::nTotStripsPerAPV
constexpr static uint16_t nTotStripsPerAPV
Definition: SiStripAPVRestorer.h:59
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
SiStripAPVRestorer::fraction_
double fraction_
Definition: SiStripAPVRestorer.h:112
SiStripAPVRestorer::derivativeFollowerRestore
void derivativeFollowerRestore(uint16_t apvN, uint16_t firstAPV, digivector_t &digis)
Definition: SiStripAPVRestorer.cc:1027
SiStripNoises::Range
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:47
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
SiStripAPVRestorer::baselineFollowerInspect
uint16_t baselineFollowerInspect(uint16_t firstAPV, const digivector_t &digis)
Definition: SiStripAPVRestorer.cc:188
SiStripQualityRcd
Definition: SiStripDependentRecords.h:56
SiStripAPVRestorer::inspectAndRestore
uint16_t inspectAndRestore(uint32_t detId, uint16_t firstAPV, const digivector_t &rawDigisPedSubtracted, digivector_t &processedRawDigi, const medians_t &vmedians)
Definition: SiStripAPVRestorer.cc:80
SiStripAPVRestorer::createCMMapCMstored
void createCMMapCMstored(const edm::DetSetVector< SiStripProcessedRawDigi > &input)
Definition: SiStripAPVRestorer.cc:1012
HcalDetIdTransform::transform
unsigned transform(const HcalDetId &id, unsigned transformCode)
Definition: HcalDetIdTransform.cc:7
SiStripAPVRestorer::flatRestore
void flatRestore(uint16_t apvN, uint16_t firstAPV, digivector_t &digis)
Definition: SiStripAPVRestorer.cc:387
SiStripAPVRestorer::meanCM_
int32_t meanCM_
Definition: SiStripAPVRestorer.h:110
SiStripQuality::IsStripBad
bool IsStripBad(const uint32_t &detid, const short &strip) const
Definition: SiStripQuality.cc:628
SiStripAPVRestorer::nSigmaNoiseDerTh_
uint32_t nSigmaNoiseDerTh_
Definition: SiStripAPVRestorer.h:116
SiStripAPVRestorer::quality_cache_id
uint32_t quality_cache_id
Definition: SiStripAPVRestorer.h:89
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
SiStripAPVRestorer::digimap_t
std::map< uint16_t, digi_t > digimap_t
Definition: SiStripAPVRestorer.h:36
SiStripPedestals::getRange
const Range getRange(const uint32_t &detID) const
Definition: SiStripPedestals.cc:28
SiStripAPVRestorer::nullInspect
uint16_t nullInspect(uint16_t firstAPV, const digivector_t &digis)
Definition: SiStripAPVRestorer.cc:303
edm::LogWarning
Definition: MessageLogger.h:141
SiStripAPVRestorer::hybridEmulationInspect
uint16_t hybridEmulationInspect(uint16_t firstAPV, const digivector_t &digis)
Definition: SiStripAPVRestorer.cc:327
SiStripAPVRestorer::apvFlags_
std::vector< std::string > apvFlags_
Definition: SiStripAPVRestorer.h:95
edm::ParameterSet
Definition: ParameterSet.h:36
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
SiStripAPVRestorer::qualityHandle
edm::ESHandle< SiStripQuality > qualityHandle
Definition: SiStripAPVRestorer.h:86
beam_dqm_sourceclient-live_cfg.cerr
cerr
Definition: beam_dqm_sourceclient-live_cfg.py:17
hgcalDigitizer_cfi.noise
noise
Definition: hgcalDigitizer_cfi.py:150
SiStripAPVRestorer::last_gradient_
int last_gradient_
Definition: SiStripAPVRestorer.h:130
pfParticleNetPreprocessParams_cfi.median
median
Definition: pfParticleNetPreprocessParams_cfi.py:16
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
SiStripAPVRestorer::apvFlagsBool_
std::vector< bool > apvFlagsBool_
Definition: SiStripAPVRestorer.h:96
SiStripAPVRestorer::filteredBaselineMax_
double filteredBaselineMax_
Definition: SiStripAPVRestorer.h:127
SiStripAPVRestorer::abnormalBaselineInspect
uint16_t abnormalBaselineInspect(uint16_t firstAPV, const digivector_t &digis)
Definition: SiStripAPVRestorer.cc:268
SiStripAPVRestorer::distortionThreshold_
uint32_t distortionThreshold_
Definition: SiStripAPVRestorer.h:119
edm::EventSetup
Definition: EventSetup.h:57
SiStripAPVRestorer::digivector_t
std::vector< digi_t > digivector_t
Definition: SiStripAPVRestorer.h:35
SiStripAPVRestorer::baselineCleaner
void baselineCleaner(const digivector_t &adcs, digimap_t &smoothedpoints, uint16_t apvN)
Definition: SiStripAPVRestorer.cc:495
SiStripAPVRestorer::size_window_
int size_window_
Definition: SiStripAPVRestorer.h:131
SiStripAPVRestorer::hitStripThreshold_
uint32_t hitStripThreshold_
Definition: SiStripAPVRestorer.h:124
get
#define get
cms::cuda::for
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
Definition: HistoContainer.h:27
SiStripAPVRestorer::restore
void restore(uint16_t firstAPV, digivector_t &digis)
Definition: SiStripAPVRestorer.cc:137
SiStripAPVRestorer::forceRestoreInspect
uint16_t forceRestoreInspect(uint16_t firstAPV, const digivector_t &digis)
Definition: SiStripAPVRestorer.cc:222
SiStripQuality::IsApvBad
bool IsApvBad(const uint32_t &detid, const short &apvNb) const
Definition: SiStripQuality.cc:620
SiStripAPVRestorer::flatRegionsFinder
bool flatRegionsFinder(const digivector_t &adcs, digimap_t &smoothedpoints, uint16_t apvN)
Definition: SiStripAPVRestorer.cc:399
SiStripAPVRestorer::cleaner_MonotonyChecker
void cleaner_MonotonyChecker(digimap_t &smoothedpoints)
Definition: SiStripAPVRestorer.cc:514
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
SiStripAPVRestorer::createCMMapRealPed
void createCMMapRealPed(const edm::DetSetVector< SiStripRawDigi > &input)
Definition: SiStripAPVRestorer.cc:993
SiStripAPVRestorer::consecThreshold_
uint32_t consecThreshold_
Definition: SiStripAPVRestorer.h:117
SiStripAPVRestorer::useRealMeanCM_
bool useRealMeanCM_
Definition: SiStripAPVRestorer.h:109
relativeConstraints.empty
bool empty
Definition: relativeConstraints.py:46
Exception
Definition: hltDiff.cc:246
SiStripAPVRestorer::restoreThreshold_
double restoreThreshold_
Definition: SiStripAPVRestorer.h:114
SiStripAPVRestorer::SiStripAPVRestorer
SiStripAPVRestorer(const edm::ParameterSet &conf)
Definition: SiStripAPVRestorer.cc:9
SiStripAPVRestorer::slopeX_
int32_t slopeX_
Definition: SiStripAPVRestorer.h:122
SiStripAPVRestorer::hybridFormatInspect
uint16_t hybridFormatInspect(uint16_t firstAPV, const digivector_t &digis)
Definition: SiStripAPVRestorer.cc:160
SiStripAPVRestorer::inspect
uint16_t inspect(uint32_t detId, uint16_t firstAPV, const digivector_t &digis, const medians_t &vmedians)
Definition: SiStripAPVRestorer.cc:90
SiStripAPVRestorer::meanCMmap_
CMMap meanCMmap_
Definition: SiStripAPVRestorer.h:92
SiStripAPVRestorer::apvFlagsBoolOverride_
std::vector< bool > apvFlagsBoolOverride_
Definition: SiStripAPVRestorer.h:97
SiStripProcessedRawDigi
A signed Digi for the silicon strip detector, containing only adc information, and suitable for stori...
Definition: SiStripProcessedRawDigi.h:16
SiStripAPVRestorer::detId_
uint32_t detId_
Definition: SiStripAPVRestorer.h:94
ztail.d
d
Definition: ztail.py:151
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
SiStripAPVRestorer::noiseHandle
edm::ESHandle< SiStripNoises > noiseHandle
Definition: SiStripAPVRestorer.h:87
SiStripPedestals::getPed
float getPed(const uint16_t &strip, const Range &range) const
Definition: SiStripPedestals.cc:51
ParameterSet.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
SiStripAPVRestorer::medians_t
std::vector< std::pair< short, float > > medians_t
Definition: SiStripAPVRestorer.h:37
SiStripAPVRestorer::baselineFollowerRestore
void baselineFollowerRestore(uint16_t apvN, uint16_t firstAPV, float median, digivector_t &digis)
Definition: SiStripAPVRestorer.cc:352
SiStripAPVRestorer::baselineAndSaturationInspect
uint16_t baselineAndSaturationInspect(uint16_t firstAPV, const digivector_t &digis)
Definition: SiStripAPVRestorer.cc:233
begin
#define begin
Definition: vmac.h:32
SiStripAPVRestorer::badAPVs_
std::vector< bool > badAPVs_
Definition: SiStripAPVRestorer.h:99
SiStripAPVRestorer.h
SiStripAPVRestorer::noise_cache_id
uint32_t noise_cache_id
Definition: SiStripAPVRestorer.h:89
SiStripAPVRestorer::applyBaselineCleaner_
bool applyBaselineCleaner_
Definition: SiStripAPVRestorer.h:120
SiStripAPVRestorer::median_
std::vector< float > median_
Definition: SiStripAPVRestorer.h:98